Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorphFuncND.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2023, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11/** \class RooMomentMorphFuncND
12 \ingroup Roofit
13
14**/
15
17
18#include <RooAbsCategory.h>
19#include <RooAbsMoment.h>
20#include <RooAddPdf.h>
21#include <RooAddition.h>
22#include <RooChangeTracker.h>
23#include <RooConstVar.h>
24#include <RooCustomizer.h>
25#include <RooFormulaVar.h>
26#include <RooLinearVar.h>
27#include <RooMoment.h>
28#include <RooNumIntConfig.h>
29#include <RooRealSumFunc.h>
30#include <RooRealVar.h>
31
33
34#include <Riostream.h>
35
36#include <TMap.h>
37#include <TMath.h>
38#include <TVector.h>
39
40#include <algorithm>
41#include <map>
42
43using std::string, std::vector;
44
45
46//_____________________________________________________________________________
47RooMomentMorphFuncND::RooMomentMorphFuncND() : _cacheMgr(this, 10, true, true), _setting(RooMomentMorphFuncND::Linear), _useHorizMorph(true)
48{
50}
51
52//_____________________________________________________________________________
56 _cacheMgr(this, 10, true, true),
57 _parList("parList", "List of morph parameters", this),
58 _obsList("obsList", "List of observables", this),
59 _referenceGrid(referenceGrid),
60 _pdfList("pdfList", "List of pdfs", this),
61 _setting(setting),
62 _useHorizMorph(true)
63{
64 // morph parameters
66
67 // observables
69
71
72 // general initialization
73 initialize();
74
76}
77
78//_____________________________________________________________________________
80 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
82 _cacheMgr(this, 10, true, true),
83 _parList("parList", "List of morph parameters", this),
84 _obsList("obsList", "List of observables", this),
85 _pdfList("pdfList", "List of pdfs", this),
86 _setting(setting),
87 _useHorizMorph(true)
88{
89 // make reference grid
90 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
92
93 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
94 for (int j = 0; j < grid.numBoundaries(); ++j) {
95 if (mrefpoints[i] == grid.array()[j]) {
96 _referenceGrid.addPdf(*static_cast<Base_t *>(pdfList.at(i)), j);
97 break;
98 }
99 }
100 }
101
103
104 // morph parameters
106 parList.add(_m);
108
109 // observables
111
112 // general initialization
113 initialize();
114
116}
117
118//_____________________________________________________________________________
120 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
122 _cacheMgr(this, 10, true, true),
123 _parList("parList", "List of morph parameters", this),
124 _obsList("obsList", "List of observables", this),
125 _pdfList("pdfList", "List of pdfs", this),
126 _setting(setting),
127 _useHorizMorph(true)
128{
129 // make reference grid
131 Int_t i = 0;
132 for (auto *mref : mrefList) {
133 if (!dynamic_cast<RooAbsReal *>(mref)) {
134 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
135 << " is not of type RooAbsReal" << std::endl;
136 throw string("RooMomentMorphFuncND::ctor() ERROR mref is not of type RooAbsReal");
137 }
138 if (!dynamic_cast<RooConstVar *>(mref)) {
139 coutW(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") WARNING mref point " << i
140 << " is not a constant, taking a snapshot of its value" << std::endl;
141 }
142 mrefpoints[i] = static_cast<RooAbsReal *>(mref)->getVal();
143 i++;
144 }
145
146 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
148 for (i = 0; i < mrefpoints.GetNrows(); ++i) {
149 for (int j = 0; j < grid.numBoundaries(); ++j) {
150 if (mrefpoints[i] == grid.array()[j]) {
151 _referenceGrid.addPdf(static_cast<Base_t &>(pdfList[i]), j);
152 break;
153 }
154 }
155 }
156
158
159 // morph parameters
161 parList.add(_m);
163
164 // observables
166
167 // general initialization
168 initialize();
169
171}
172
173//_____________________________________________________________________________
176 _cacheMgr(other._cacheMgr, this),
177 _parList("parList", this, other._parList),
178 _obsList("obsList", this, other._obsList),
179 _referenceGrid(other._referenceGrid),
180 _pdfList("pdfList", this, other._pdfList),
181 _setting(other._setting),
182 _useHorizMorph(other._useHorizMorph),
183 _isPdfMode{other._isPdfMode}
184{
185 // general initialization
186 initialize();
187
189}
190
191//_____________________________________________________________________________
196
197//_____________________________________________________________________________
199{
200 for (vector<RooAbsBinning *>::iterator itr = _referenceGrid._grid.begin(); itr != _referenceGrid._grid.end();
201 ++itr) {
202 _referenceGrid._nnuis.push_back((*itr)->numBins() + 1);
203 }
204
205 int nPar = _parList.size();
206 int nDim = _referenceGrid._grid.size();
208 int nRef = _referenceGrid._nref.size();
209 int depth = std::pow(2, nPar);
210
211 if (nPar != nDim) {
212 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPar != nDim"
213 << ": " << nPar << " !=" << nDim << std::endl;
214 assert(0);
215 }
216
217 if (nPdf != nRef) {
218 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPdf != nRef"
219 << ": " << nPdf << " !=" << nRef << std::endl;
220 assert(0);
221 }
222
223 // Transformation matrix for NonLinear settings
224 _M = std::make_unique<TMatrixD>(nPdf, nPdf);
225 _MSqr = std::make_unique<TMatrixD>(depth, depth);
227 TMatrixD M(nPdf, nPdf);
228
230 for (int k = 0; k < nPdf; ++k) {
232 for (int idim = 0; idim < nPar; idim++) {
233 double delta = _referenceGrid._nref[k][idim] - _referenceGrid._nref[0][idim];
234 dm2.push_back(delta);
235 }
236 dm[k] = dm2;
237 }
238
240 for (int idim = 0; idim < nPar; idim++) {
242 xtmp.reserve(_referenceGrid._nnuis[idim]);
243 for (int ix = 0; ix < _referenceGrid._nnuis[idim]; ix++) {
244 xtmp.push_back(ix);
245 }
246 powers.push_back(xtmp);
247 }
248
249 vector<vector<int>> output;
251 int nCombs = output.size();
252
253 for (int k = 0; k < nPdf; ++k) {
254 int nperm = 0;
255 for (int i = 0; i < nCombs; i++) {
256 double tmpDm = 1.0;
257 for (int ix = 0; ix < nPar; ix++) {
258 double delta = dm[k][ix];
259 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
260 }
261 M(k, nperm) = tmpDm;
262 nperm++;
263 }
264 }
265
266 // M.Print();
267 (*_M) = M.Invert();
268 }
269
270 // Resize transformation vectors
271 _squareVec.resize(std::pow(2, nPar));
272 _squareIdx.resize(std::pow(2, nPar));
273}
274
275//_____________________________________________________________________________
277 : _pdfList(other._pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
278{
279 for (unsigned int i = 0; i < other._grid.size(); i++) {
280 _grid.push_back(other._grid[i]->clone());
281 }
282}
283
284//_____________________________________________________________________________
286{
287 for (RooAbsBinning *binning : _grid) {
288 delete binning;
289 }
290}
291
292//_____________________________________________________________________________
294{
297 thisBoundaries.push_back(bin_x);
298 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
299 _pdfList.add(pdf);
300 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
301 _nref.push_back(thisBoundaryCoordinates);
302}
303
304//_____________________________________________________________________________
306{
309 thisBoundaries.push_back(bin_x);
310 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
311 thisBoundaries.push_back(bin_y);
312 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
313 _pdfList.add(pdf);
314 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
315 _nref.push_back(thisBoundaryCoordinates);
316}
317
318//_____________________________________________________________________________
320{
323 thisBoundaries.push_back(bin_x);
324 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
325 thisBoundaries.push_back(bin_y);
326 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
327 thisBoundaries.push_back(bin_z);
328 thisBoundaryCoordinates.push_back(_grid[2]->array()[bin_z]);
329 _pdfList.add(pdf);
330 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
331 _nref.push_back(thisBoundaryCoordinates);
332}
333
334//_____________________________________________________________________________
336{
338 int nBins = bins.size();
339 thisBoundaryCoordinates.reserve(nBins);
340 for (int i = 0; i < nBins; i++) {
341 thisBoundaryCoordinates.push_back(_grid[i]->array()[bins[i]]);
342 }
343 _pdfList.add(pdf);
344 _pdfMap[bins] = _pdfList.size() - 1;
345 _nref.push_back(thisBoundaryCoordinates);
346}
347
348//_____________________________________________________________________________
349std::unique_ptr<RooAbsArg>
351{
352 // Build (or fetch) the cache that holds the morph's internal compute graph:
353 // moment integrals, slope/offset formulas, RooLinearVar transforms, the per-pdf
354 // RooHistPdf clones, and the final RooAddPdf/RooRealSumFunc sum.
355 CacheElem *cache = getCache(&normSet);
356
357 // Make sure fractions hold sensible initial values (the replacement nodes
358 // below will keep them in sync going forward).
359 cache->calculateFractions(*this, false);
360
361 // The cache's subtree carries ORIGNAME: attributes left over from the
362 // RooCustomizer that built the per-pdf transformed RooHistPdf clones (it
363 // tagged each transVar with ORIGNAME:<obs>). RooCustomizer::build calls
364 // redirectServers with nameChange=true and will throw if it sees several
365 // candidates with the same ORIGNAME:* attribute, which is exactly what
366 // happens here once we ask it to clone the subtree. Strip those stale
367 // markers before cloning.
368 {
370 cache->_sum->branchNodeServerList(&branches);
371 for (auto *b : branches) {
372 std::vector<std::string> toRemove;
373 for (auto const &attr : b->attributes()) {
374 if (attr.rfind("ORIGNAME:", 0) == 0)
375 toRemove.push_back(attr);
376 }
377 for (auto const &attr : toRemove)
378 b->setAttribute(attr.c_str(), false);
379 }
380 }
381
382 // Replace each of the imperatively-updated fraction RooRealVars with a
383 // RooMomentMorphFraction node. This puts the fraction-recomputation inside
384 // the Evaluator's compute graph, so the moment integrals (which the
385 // fractions depend on transitively via the slope/offset formulas) become
386 // sibling nodes that the Evaluator caches once per minimization step.
388 const int nFrac = cache->_frac.size();
389 for (int i = 0; i < nFrac; ++i) {
390 auto frac = static_cast<RooRealVar *>(cache->_frac.at(i));
391 std::string newName = std::string{frac->GetName()} + "_compiled";
392 newFractions.addOwned(
393 std::make_unique<RooFit::Detail::RooMomentMorphFraction>(newName.c_str(), frac->GetTitle(), *this, i));
394 }
395
397 RooCustomizer cust(*cache->_sum, "compiled");
398 cust.setCloneBranchSet(clonedBranches);
399 for (int i = 0; i < nFrac; ++i) {
400 cust.replaceArg(*cache->_frac.at(i), newFractions[i]);
401 }
402
403 // RooCustomizer::build() already transfers ownership of the cloned branches
404 // (everything in `clonedBranches` except the returned top node) to the new
405 // top node's owned-components list, so we only have to attach the
406 // newly-created fraction nodes here.
407 std::unique_ptr<RooAbsReal> newSum{static_cast<RooAbsReal *>(cust.build())};
408 newSum->addOwnedComponents(std::move(newFractions));
409
410 if (_isPdfMode) {
411 // In pdf mode, _sum is a RooAddPdf whose value should be the normalized
412 // morph density. We must let the inner RooHistPdf clones go through
413 // their own compileForNormSet so they get wrapped in a RooNormalizedPdf.
415 newSum->branchNodeServerList(&branches);
416 for (auto *b : branches) {
417 if (dynamic_cast<RooFit::Detail::RooMomentMorphFraction *>(b))
418 ctx.markAsCompiled(*b);
419 }
421 } else {
422 // Non-pdf mode: _sum is a RooRealSumFunc and the legacy path returns
423 // the raw weighted sum of the per-pdf bin values, with no per-component
424 // normalization. So cloned subtree is already in its final form.
427 }
428
429 return newSum;
430}
431
432namespace RooFit {
433namespace Detail {
434
435RooMomentMorphFraction::RooMomentMorphFraction(const char *name, const char *title, RooMomentMorphFuncND const &parent,
436 int index)
437 : RooAbsReal(name, title), _parList("parList", "parList", this), _parent(&parent), _index(index)
438{
439 _parList.add(parent._parList);
440}
441
443 : RooAbsReal(other, name), _parList("parList", this, other._parList), _parent(other._parent), _index(other._index)
444{
445}
446
448{
449 auto *cache = _parent->getCache(nullptr);
450 if (cache->_tracker->hasChanged(true)) {
451 cache->calculateFractions(*_parent, false);
452 }
453 return cache->frac(_index)->getVal();
454}
455
456} // namespace Detail
457} // namespace RooFit
458
459//_____________________________________________________________________________
461{
462 auto cache = static_cast<CacheElem *>(_cacheMgr.getObj(nullptr, static_cast<RooArgSet const*>(nullptr)));
463 if (cache) {
464 return cache;
465 }
466
467 int nObs = _obsList.size();
469
470 RooAbsReal *null = nullptr;
471 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
472 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
473 vector<RooAbsReal *> myrms(nObs, null);
474 vector<RooAbsReal *> mypos(nObs, null);
475 vector<RooAbsReal *> slope(nPdf * nObs, null);
476 vector<RooAbsReal *> offsets(nPdf * nObs, null);
477 vector<RooAbsReal *> transVar(nPdf * nObs, null);
479
482
483 // fraction parameters
484 RooArgList coefList("coefList"); // fractions multiplied with input pdfs
485 RooArgList coefList2("coefList2"); // fractions multiplied with mean position of observable contribution
486 RooArgList coefList3("coefList3"); // fractions multiplied with rms position of observable contribution
487
488 for (int i = 0; i < 3 * nPdf; ++i) {
489 string fracName = Form("frac_%d", i);
490 double initval = _isPdfMode ? 1.0 : 0.0;
491 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), initval); // to be set later
492
493 fracl.add(*frac);
494 if (i < nPdf) {
495 coefList.add(*static_cast<RooRealVar *>(fracl.at(i)));
496 } else if (i < 2 * nPdf) {
497 coefList2.add(*static_cast<RooRealVar *>(fracl.at(i)));
498 } else {
499 coefList3.add(*static_cast<RooRealVar *>(fracl.at(i)));
500 }
501 ownedComps.add(*static_cast<RooRealVar *>(fracl.at(i)));
502 }
503
504 std::unique_ptr<RooAbsReal> theSum;
505 string sumName = Form("%s_sum", GetName());
506
508 if (_useHorizMorph) {
509 // mean and sigma
511 for (int i = 0; i < nPdf; ++i) {
512 for (int j = 0; j < nObs; ++j) {
513 RooAbsMoment *mom = nObs == 1 ? (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)))
514 : (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)), obsList);
515
516 mom->setLocalNoDirtyInhibit(true);
517 mom->mean()->setLocalNoDirtyInhibit(true);
518
519 sigmarv[sij(i, j)] = mom;
520 meanrv[sij(i, j)] = mom->mean();
521
522 ownedComps.add(*sigmarv[sij(i, j)]);
523 }
524 }
525
526 // slope and offset (to be set later, depend on nuisance parameters)
527 for (int j = 0; j < nObs; ++j) {
528 RooArgList meanList("meanList");
529 RooArgList rmsList("rmsList");
530 for (int i = 0; i < nPdf; ++i) {
531 meanList.add(*meanrv[sij(i, j)]);
532 rmsList.add(*sigmarv[sij(i, j)]);
533 }
534 string myrmsName = Form("%s_rms_%d", GetName(), j);
535 string myposName = Form("%s_pos_%d", GetName(), j);
536 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
537 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
538 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
539 }
540
541 // construction of unit pdfs
542
543 Int_t i = 0;
544 for (auto const *pdf : static_range_cast<Base_t *>(_pdfList)) {
545
546 string pdfName = Form("pdf_%d", i);
547 RooCustomizer cust(*pdf, pdfName.c_str());
548
549 Int_t j = 0;
550 for (auto *var : static_range_cast<RooRealVar *>(obsList)) {
551 // slope and offset formulas
552 string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
553 string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
554
555 slope[sij(i, j)] =
556 new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[sij(i, j)], *myrms[j]));
557 offsets[sij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
558 RooArgList(*meanrv[sij(i, j)], *mypos[j], *slope[sij(i, j)]));
559 ownedComps.add(RooArgSet(*slope[sij(i, j)], *offsets[sij(i, j)]));
560
561 // linear transformations, so pdf can be renormalized easily
562 string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
563 transVar[sij(i, j)] = new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[sij(i, j)],
564 *offsets[sij(i, j)]);
565
566 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
567 // This will prevent the likelihood optimizers from erroneously declaring terms constant
568 transVar[sij(i, j)]->addServerList((RooAbsCollection &)_parList);
569
570 ownedComps.add(*transVar[sij(i, j)]);
571 cust.replaceArg(*var, *transVar[sij(i, j)]);
572 ++j;
573 }
574 transPdf[i] = static_cast<Base_t *>(cust.build());
575 transPdfList.add(*transPdf[i]);
576 ownedComps.add(*transPdf[i]);
577 ++i;
578 }
579 }
580
581 // sum pdf
582 RooArgList const &pdfList = _useHorizMorph ? transPdfList : static_cast<RooArgList const &>(_pdfList);
583 if (_isPdfMode) {
584 theSum = std::make_unique<RooAddPdf>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
585 } else {
586 theSum = std::make_unique<RooRealSumFunc>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
587 }
588
589 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
590 // This will prevent the likelihood optimizers from erroneously declaring terms constant
591 theSum->addServerList((RooAbsCollection &)_parList);
592 theSum->addOwnedComponents(ownedComps);
593
594 // change tracker for fraction parameters
595 std::string trackerName = std::string(GetName()) + "_frac_tracker";
596
597 // Store it in the cache
598 cache = new CacheElem(std::move(theSum),
599 std::make_unique<RooChangeTracker>(trackerName.c_str(), trackerName.c_str(), _parList, true),
600 fracl);
601 _cacheMgr.setObj(nullptr, nullptr, cache, nullptr);
602
603 return cache;
604}
605
606RooMomentMorphFuncND::CacheElem::CacheElem(std::unique_ptr<RooAbsReal> &&sumFunc,
607 std::unique_ptr<RooChangeTracker> &&tracker, const RooArgList &flist)
608 : _sum(std::move(sumFunc)), _tracker(std::move(tracker))
609{
610 _frac.add(flist);
611}
612
613//_____________________________________________________________________________
618
619//_____________________________________________________________________________
621
622//_____________________________________________________________________________
624{
625 // Special version of getValV() overrides Base_t::getValV() to save value of current normalization set
626 _curNormSet = set ? const_cast<RooArgSet *>(set) : const_cast<RooArgSet *>(static_cast<RooArgSet const*>(&_obsList));
627 return Base_t::getValV(set);
628}
629
630//_____________________________________________________________________________
632{
633 CacheElem *cache = getCache(nset ? nset : _curNormSet);
634
635 if (cache->_tracker->hasChanged(true)) {
636 cache->calculateFractions(*this, false); // verbose turned off
637 }
638 return cache->_sum.get();
639}
640
641//_____________________________________________________________________________
643{
645
646 if (cache->_tracker->hasChanged(true)) {
647 cache->calculateFractions(*this, false); // verbose turned off
648 }
649
650 double ret = cache->_sum->getVal(_obsList.nset());
651
652 return ret;
653}
654
655//_____________________________________________________________________________
657{
658 return static_cast<RooRealVar *>(_frac.at(i));
659}
660
661//_____________________________________________________________________________
663{
664 return static_cast<RooRealVar *>(_frac.at(i));
665}
666
667//_____________________________________________________________________________
669{
670 int nPdf = self._pdfList.size();
671 int nPar = self._parList.size();
672
673 double fracLinear(1.);
674 double fracNonLinear(1.);
675
676 if (self._setting == NonLinear || self._setting == NonLinearLinFractions || self._setting == NonLinearPosFractions) {
677 // Calculate the delta vector
679 for (int idim = 0; idim < nPar; idim++) {
680 double delta = (static_cast<RooRealVar *>(self._parList.at(idim)))->getVal() - self._referenceGrid._nref[0][idim];
681 dm2.push_back(delta);
682 }
683
685 for (int idim = 0; idim < nPar; idim++) {
687 xtmp.reserve(self._referenceGrid._nnuis[idim]);
688 for (int ix = 0; ix < self._referenceGrid._nnuis[idim]; ix++) {
689 xtmp.push_back(ix);
690 }
691 powers.push_back(xtmp);
692 }
693
694 vector<vector<int>> output;
696 int nCombs = output.size();
697
699
700 int nperm = 0;
701 for (int i = 0; i < nCombs; i++) {
702 double tmpDm = 1.0;
703 for (int ix = 0; ix < nPar; ix++) {
704 double delta = dm2[ix];
705 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
706 }
708 nperm++;
709 }
710
711 double sumposfrac = 0.0;
712 for (int i = 0; i < nPdf; ++i) {
713 double ffrac = 0.0;
714
715 for (int j = 0; j < nPdf; ++j) {
716 ffrac += (*self._M)(j, i) * deltavec[j] * fracNonLinear;
717 }
718
719 if (ffrac >= 0) {
720 sumposfrac += ffrac;
721 }
722
723 // fractions for pdf
724 if (self._setting != NonLinearLinFractions) {
725 const_cast<RooRealVar *>(frac(i))->setVal(ffrac);
726 }
727
728 // fractions for rms and mean
729 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(ffrac); // need to add up
730 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(ffrac); // need to add up
731
732 if (verbose) {
733 std::cout << "NonLinear fraction " << ffrac << std::endl;
734 frac(i)->Print();
735 frac(nPdf + i)->Print();
736 frac(2 * nPdf + i)->Print();
737 }
738 }
739
740 if (self._setting == NonLinearPosFractions) {
741 for (int i = 0; i < nPdf; ++i) {
742 if (frac(i)->getVal() < 0)
743 const_cast<RooRealVar *>(frac(i))->setVal(0.);
744 const_cast<RooRealVar *>(frac(i))->setVal(frac(i)->getVal() / sumposfrac);
745 }
746 }
747 }
748
749 if (self._setting == Linear || self._setting == NonLinearLinFractions) {
750 // zero all fractions
751 // for (int i = 0; i < 3*nPdf; ++i) {
752 for (int i = 0; i < nPdf; ++i) {
753 double initval = 0;
754 const_cast<RooRealVar *>(frac(i))->setVal(initval);
755 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(initval);
756 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(initval);
757 }
758
759 std::vector<double> mtmp;
760
761 // loop over parList
762 for (auto *m : static_range_cast<RooRealVar *>(self._parList)) {
763 mtmp.push_back(m->getVal());
764 }
765
766 self.findShape(mtmp); // this sets _squareVec and _squareIdx quantities
767
768 int depth = std::pow(2, nPar);
770
771 int nperm = 0;
772
774 xtmp.reserve(nPar);
775 for (int ix = 0; ix < nPar; ix++) {
776 xtmp.push_back(ix);
777 }
778
779 for (int iperm = 1; iperm <= nPar; ++iperm) {
780 do {
781 double dtmp = mtmp[xtmp[0]] - self._squareVec[0][xtmp[0]];
782 for (int itmp = 1; itmp < iperm; ++itmp) {
783 dtmp *= mtmp[xtmp[itmp]] - self._squareVec[0][xtmp[itmp]];
784 }
785 deltavec[nperm + 1] = dtmp;
786 nperm++;
788 }
789
790 double origFrac1(0.);
791 double origFrac2(0.);
792 for (int i = 0; i < depth; ++i) {
793 double ffrac = 0.;
794 for (int j = 0; j < depth; ++j) {
795 ffrac += (*self._MSqr)(j, i) * deltavec[j] * fracLinear;
796 }
797
798 // set fractions for pdf
799 origFrac1 = frac(self._squareIdx[i])->getVal(); // already set in case of smoothlinear
800 const_cast<RooRealVar *>(frac(self._squareIdx[i]))->setVal(origFrac1 + ffrac); // need to add up
801
802 // set fractions for rms and mean
803 if (self._setting != NonLinearLinFractions) {
804 origFrac2 =
805 frac(nPdf + self._squareIdx[i])->getVal(); // already set in case of smoothlinear
806 const_cast<RooRealVar *>(frac(nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
807 const_cast<RooRealVar *>(frac(2 * nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
808 }
809
810 if (verbose) {
811 std::cout << "Linear fraction " << ffrac << std::endl;
812 frac(self._squareIdx[i])->Print();
813 frac(nPdf + self._squareIdx[i])->Print();
814 frac(2 * nPdf + self._squareIdx[i])->Print();
815 }
816 }
817 }
818}
819
820//_____________________________________________________________________________
822{
823 int nPar = _parList.size();
824 int nRef = _referenceGrid._nref.size();
825
826 // Find hypercube enclosing the location to morph to
827 // bool isEnclosed = true;
828 // for (int i = 0; i < nPar; i++) {
829 // if (x[i] < _referenceGrid._grid[i]->lowBound())
830 // isEnclosed = false;
831 // if (x[i] > _referenceGrid._grid[i]->highBound())
832 // isEnclosed = false;
833 // }
834
835 // std::cout << "isEnclosed = " << isEnclosed << std::endl;
836
837 int depth = std::pow(2, nPar);
838
839 vector<vector<double>> boundaries(nPar);
840 for (int idim = 0; idim < nPar; idim++) {
841 int bin = _referenceGrid._grid[idim]->binNumber(x[idim]);
842 double lo = _referenceGrid._grid[idim]->binLow(bin);
843 double hi = _referenceGrid._grid[idim]->binHigh(bin);
844 boundaries[idim].push_back(lo);
845 boundaries[idim].push_back(hi);
846 }
847
848 vector<vector<double>> output;
849 RooFit::Detail::cartesianProduct(output, boundaries);
850 _squareVec = output;
851
852 for (int isq = 0; isq < depth; isq++) {
853 for (int iref = 0; iref < nRef; iref++) {
856 break;
857 }
858 }
859 }
860
861 // std::cout << std::endl;
862
863 // for (int isq = 0; isq < _squareVec.size(); isq++) {
864 // std::cout << _squareIdx[isq];
865 // std::cout << " (";
866 // for (int isqq = 0; isqq < _squareVec[isq].size(); isqq++) {
867 // std::cout << _squareVec[isq][isqq] << ((isqq<_squareVec[isq].size()-1)?",":"");
868 // }
869 // std::cout << ") ";
870 // }
871
872 // construct transformation matrix for linear extrapolation
873 TMatrixD M(depth, depth);
874
876 xtmp.reserve(nPar);
877 for (int ix = 0; ix < nPar; ix++) {
878 xtmp.push_back(ix);
879 }
880
881 for (int k = 0; k < depth; ++k) {
882 M(k, 0) = 1.0;
883
884 int nperm = 0;
886
887 for (int iperm = 1; iperm <= nPar; ++iperm) {
888 do {
889 double dtmp = _squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
890 for (int itmp = 1; itmp < iperm; ++itmp) {
892 }
893 M(k, nperm + 1) = dtmp;
894 nperm++;
896 }
897 }
898
899 // M.Print();
900 (*_MSqr) = M.Invert();
901}
902
903//_____________________________________________________________________________
905{
906 if (allVars.size() == 1) {
907 RooAbsReal *temp = const_cast<RooMomentMorphFuncND *>(this);
908 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator");
909 int nbins = (static_cast<RooRealVar *>(allVars.first()))->numBins();
910 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins", nbins);
911 return true;
912 } else {
913 std::cout << "Currently BinIntegrator only knows how to deal with 1-d " << std::endl;
914 return false;
915 }
916 return false;
917}
#define b(i)
Definition RSha256.hxx:100
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t attr
char name[80]
Definition TGX11.cxx:148
#define hi
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2496
const_iterator begin() const
const_iterator end() const
Abstract base class for RooRealVar binning definitions.
Abstract container object that can hold multiple RooAbsArg objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * first() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
virtual double getValV(const RooArgSet *normalisationSet=nullptr) const
Return value of object.
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition RooAddition.h:27
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Implements a RooAbsBinning in terms of an array of boundary values, posing no constraints on the choi...
Definition RooBinning.h:27
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
bool setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Represents a constant real-valued object.
Definition RooConstVar.h:23
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void markAsCompiled(RooAbsArg &arg) const
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
void markSubtreeAsCompiled(RooAbsArg &arg) const
Mark arg and every branch node reachable through its server tree as already compiled.
Helper compute-graph node that exposes one of the morph mixing fractions to the RooFit::Evaluator.
const RooMomentMorphFuncND * _parent
! morph that owns the cache (not owned)
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooLinearVar is the most general form of a derived real-valued object that can be used by RooRealInte...
void calculateFractions(const RooMomentMorphFuncND &self, bool verbose=true) const
std::unique_ptr< RooChangeTracker > _tracker
std::unique_ptr< RooAbsReal > _sum
RooArgList containedArgs(Action) override
CacheElem(std::unique_ptr< RooAbsReal > &&sumFunc, std::unique_ptr< RooChangeTracker > &&tracker, const RooArgList &flist)
void addBinning(const RooAbsBinning &binning)
std::vector< RooAbsBinning * > _grid
std::vector< std::vector< double > > _nref
void addPdf(const RooAbsReal &func, int bin_x)
RooObjCacheManager _cacheMgr
! Transient cache manager
std::unique_ptr< TMatrixD > _MSqr
void findShape(const std::vector< double > &x) const
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsReal * sumFunc(const RooArgSet *nset)
CacheElem * getCache(const RooArgSet *nset) const
std::unique_ptr< TMatrixD > _M
RooArgSet * _curNormSet
! Transient cache manager
bool setBinIntegrator(RooArgSet &allVars)
double getValV(const RooArgSet *set=nullptr) const override
Return value of object.
std::vector< std::vector< double > > _squareVec
int sij(const int &i, const int &j) const
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
std::vector< int > _squareIdx
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooCategory & method1D()
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Double_t x[n]
Definition legend1.C:17
void cartesianProduct(std::vector< std::vector< T > > &out, std::vector< std::vector< T > > &in)
Definition Algorithms.h:22
bool nextCombination(const Iterator first, Iterator k, const Iterator last)
Definition Algorithms.h:64
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:72
TMarker m
Definition textangle.C:8