Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorphND.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * This code was autogenerated by RooClassFactory *
5 *****************************************************************************/
6
7#include "Riostream.h"
8
9#include "RooMomentMorphND.h"
10#include "RooAbsCategory.h"
11#include "RooRealIntegral.h"
12#include "RooRealConstant.h"
13#include "RooRealVar.h"
14#include "RooFormulaVar.h"
15#include "RooCustomizer.h"
16#include "RooAddPdf.h"
17#include "RooAddition.h"
18#include "RooAbsMoment.h"
19#include "RooMoment.h"
20#include "RooLinearVar.h"
21#include "RooChangeTracker.h"
22#include "RooNumIntConfig.h"
23#include "RooHistPdf.h"
24
25#include "TMath.h"
26#include "TVector.h"
27#include "TMap.h"
28
29#include <map>
30#include <algorithm>
31
32using namespace std;
33
35
36//_____________________________________________________________________________
38 : _cacheMgr(this, 10, true, true), _curNormSet(0), _M(0), _MSqr(0), _setting(RooMomentMorphND::Linear), _useHorizMorph(true)
39{
40 _parItr = _parList.createIterator();
41 _obsItr = _obsList.createIterator();
42
44}
45
46//_____________________________________________________________________________
47RooMomentMorphND::RooMomentMorphND(const char *name, const char *title, const RooArgList &parList,
48 const RooArgList &obsList, const Grid &referenceGrid, const Setting &setting)
49 : RooAbsPdf(name, title), _cacheMgr(this, 10, kTRUE, kTRUE), _parList("parList", "List of morph parameters", this),
50 _obsList("obsList", "List of observables", this), _referenceGrid(referenceGrid),
51 _pdfList("pdfList", "List of pdfs", this), _setting(setting), _useHorizMorph(true)
52{
53 // morph parameters
54 initializeParameters(parList);
55
56 // observables
57 initializeObservables(obsList);
58
60
61 // general initialization
62 initialize();
63
65}
66
67//_____________________________________________________________________________
68RooMomentMorphND::RooMomentMorphND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
69 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
70 : RooAbsPdf(name, title), _cacheMgr(this, 10, kTRUE, kTRUE), _parList("parList", "List of morph parameters", this),
71 _obsList("obsList", "List of observables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
72 _useHorizMorph(true)
73{
74 // make reference grid
75 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
77
78 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
79 for (int j = 0; j < grid.numBoundaries(); ++j) {
80 if (mrefpoints[i] == grid.array()[j]) {
81 _referenceGrid.addPdf(*(RooAbsPdf *)pdfList.at(i), j);
82 break;
83 }
84 }
85 }
86
88
89 // morph parameters
90 RooArgList parList;
91 parList.add(_m);
92 initializeParameters(parList);
93
94 // observables
95 initializeObservables(varList);
96
97 // general initialization
98 initialize();
99
101}
102
103//_____________________________________________________________________________
104RooMomentMorphND::RooMomentMorphND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
105 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
106 : RooAbsPdf(name, title), _cacheMgr(this, 10, kTRUE, kTRUE), _parList("parList", "List of morph parameters", this),
107 _obsList("obsList", "List of observables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
108 _useHorizMorph(true)
109{
110 // make reference grid
111 TVectorD mrefpoints(mrefList.getSize());
112 TIterator *mrefItr = mrefList.createIterator();
113 RooAbsReal *mref;
114 for (int i = 0; (mref = dynamic_cast<RooAbsReal *>(mrefItr->Next())); ++i) {
115 if (!mref) {
116 coutE(InputArguments) << "RooMomentMorphND::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
117 << " is not of type RooAbsReal" << endl;
118 throw string("RooMomentMorphND::ctor() ERROR mref is not of type RooAbsReal");
119 }
120 if (!dynamic_cast<RooConstVar *>(mref)) {
121 coutW(InputArguments) << "RooMomentMorphND::ctor(" << GetName() << ") WARNING mref point " << i
122 << " is not a constant, taking a snapshot of its value" << endl;
123 }
124 mrefpoints[i] = mref->getVal();
125 }
126 delete mrefItr;
127
128 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
130
131 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
132 for (int j = 0; j < grid.numBoundaries(); ++j) {
133 if (mrefpoints[i] == grid.array()[j]) {
134 _referenceGrid.addPdf(*(RooAbsPdf *)pdfList.at(i), j);
135 break;
136 }
137 }
138 }
139
141
142 // morph parameters
143 RooArgList parList;
144 parList.add(_m);
145 initializeParameters(parList);
146
147 // observables
148 initializeObservables(varList);
149
150 // general initialization
151 initialize();
152
154}
155
156//_____________________________________________________________________________
158 : RooAbsPdf(other, name), _cacheMgr(other._cacheMgr, this), _curNormSet(0),
159 _parList("parList", this, other._parList), _obsList("obsList", this, other._obsList),
160 _referenceGrid(other._referenceGrid), _pdfList("pdfList", this, other._pdfList), _M(0), _MSqr(0),
161 _setting(other._setting), _useHorizMorph(other._useHorizMorph)
162{
165
166 // general initialization
167 initialize();
168
170}
171
172//_____________________________________________________________________________
174{
175 if (_parItr)
176 delete _parItr;
177 if (_obsItr)
178 delete _obsItr;
179 if (_M)
180 delete _M;
181 if (_MSqr)
182 delete _MSqr;
183
185}
186
187//_____________________________________________________________________________
189{
190 TIterator *parItr = parList.createIterator();
191 RooAbsArg *par;
192 for (int i = 0; (par = (RooAbsArg *)parItr->Next()); ++i) {
193 if (!dynamic_cast<RooAbsReal *>(par)) {
194 coutE(InputArguments) << "RooMomentMorphND::ctor(" << GetName() << ") ERROR: parameter " << par->GetName()
195 << " is not of type RooAbsReal" << endl;
196 throw string("RooMomentMorphND::initializeParameters() ERROR parameter is not of type RooAbsReal");
197 }
198 _parList.add(*par);
199 }
200 delete parItr;
201
203}
204
205//_____________________________________________________________________________
207{
208 TIterator *obsItr = obsList.createIterator();
209 RooAbsArg *var;
210 for (int i = 0; (var = (RooAbsArg *)obsItr->Next()); ++i) {
211 if (!dynamic_cast<RooAbsReal *>(var)) {
212 coutE(InputArguments) << "RooMomentMorphND::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
213 << " is not of type RooAbsReal" << endl;
214 throw string("RooMomentMorphND::initializeObservables() ERROR variable is not of type RooAbsReal");
215 }
216 _obsList.add(*var);
217 }
218 delete obsItr;
219
221}
222
223//_____________________________________________________________________________
224// from http://stackoverflow.com/a/5279601
225template <typename T>
226struct Digits {
227 typename vector<T>::const_iterator begin;
228 typename vector<T>::const_iterator end;
229 typename vector<T>::const_iterator me;
230};
231
232template <typename T>
233inline void cartesian_product(vector<vector<T>> &out, vector<vector<T>> &in)
234{
235 vector<Digits<T>> vd;
236
237 for (typename vector<vector<T>>::const_iterator it = in.begin(); it != in.end(); ++it) {
238 Digits<T> d = {(*it).begin(), (*it).end(), (*it).begin()};
239 vd.push_back(d);
240 }
241
242 while (1) {
243 vector<T> result;
244 for (typename vector<Digits<T>>::const_iterator it = vd.begin(); it != vd.end(); ++it) {
245 result.push_back(*(it->me));
246 }
247 out.push_back(result);
248
249 for (typename vector<Digits<T>>::iterator it = vd.begin();;) {
250 ++(it->me);
251 if (it->me == it->end) {
252 if (it + 1 == vd.end()) {
253 return;
254 } else {
255 it->me = it->begin;
256 ++it;
257 }
258 } else {
259 break;
260 }
261 }
262 }
263}
264
265//_____________________________________________________________________________
267{
268 // TIterator* pdfItr = _referenceGrid._pdfList.createIterator() ;
269 // RooAbsPdf* pdf ;
270 // for (int i=0; (pdf = dynamic_cast<RooAbsPdf*>(pdfItr->Next())); ++i) {
271 // if (!pdf) {
272 // coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: pdf " << pdf->GetName() << " is not
273 // of type RooAbsPdf" << endl ;
274 // throw string("RooPolyMorh::ctor() ERROR pdf is not of type RooAbsPdf") ;
275 // }
276 // _pdfList.addClone(*pdf) ;
277 // }
278 // delete pdfItr ;
279
280 for (vector<RooAbsBinning *>::iterator itr = _referenceGrid._grid.begin(); itr != _referenceGrid._grid.end();
281 ++itr) {
282 _referenceGrid._nnuis.push_back((*itr)->numBins() + 1);
283 }
284
285 int nPar = _parList.getSize();
286 int nDim = _referenceGrid._grid.size();
287 int nPdf = _referenceGrid._pdfList.getSize();
288 int nRef = _referenceGrid._nref.size();
289 int depth = TMath::Power(2, nPar);
290
291 if (nPar != nDim) {
292 coutE(InputArguments) << "RooMomentMorphND::initialize(" << GetName() << ") ERROR: nPar != nDim"
293 << ": " << nPar << " !=" << nDim << endl;
294 assert(0);
295 }
296
297 if (nPdf != nRef) {
298 coutE(InputArguments) << "RooMomentMorphND::initialize(" << GetName() << ") ERROR: nPdf != nRef"
299 << ": " << nPdf << " !=" << nRef << endl;
300 assert(0);
301 }
302
303 // Transformation matrix for NonLinear settings
304 _M = new TMatrixD(nPdf, nPdf);
305 _MSqr = new TMatrixD(depth, depth);
307 TMatrixD M(nPdf, nPdf);
308
309 vector<vector<double>> dm(nPdf);
310 for (int k = 0; k < nPdf; ++k) {
311 vector<double> dm2;
312 for (int idim = 0; idim < nPar; idim++) {
313 Double_t delta = _referenceGrid._nref[k][idim] - _referenceGrid._nref[0][idim];
314 dm2.push_back(delta);
315 }
316 dm[k] = dm2;
317 }
318
319 vector<vector<int>> powers;
320 for (int idim = 0; idim < nPar; idim++) {
321 vector<int> xtmp;
322 for (int ix = 0; ix < _referenceGrid._nnuis[idim]; ix++) {
323 xtmp.push_back(ix);
324 }
325 powers.push_back(xtmp);
326 }
327
328 vector<vector<int>> output;
329 cartesian_product(output, powers);
330 int nCombs = output.size();
331
332 for (int k = 0; k < nPdf; ++k) {
333 int nperm = 0;
334 for (int i = 0; i < nCombs; i++) {
335 double tmpDm = 1.0;
336 for (int ix = 0; ix < nPar; ix++) {
337 Double_t delta = dm[k][ix];
338 tmpDm *= TMath::Power(delta, static_cast<double>(output[i][ix]));
339 }
340 M(k, nperm) = tmpDm;
341 nperm++;
342 }
343 }
344
345 // M.Print();
346 (*_M) = M.Invert();
347 }
348
349 // Resize transformation vectors
350 _squareVec.resize(TMath::Power(2, nPar));
351 _squareIdx.resize(TMath::Power(2, nPar));
352}
353
354//_____________________________________________________________________________
356 : _grid(other._grid), _pdfList(other._pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
357{
358}
359
360//_____________________________________________________________________________
362{
363}
364
365//_____________________________________________________________________________
366void RooMomentMorphND::Grid::addPdf(const RooAbsPdf &pdf, int bin_x)
367{
368 vector<int> thisBoundaries;
369 vector<double> thisBoundaryCoordinates;
370 thisBoundaries.push_back(bin_x);
371 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
372 _pdfList.add(pdf);
373 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
374 _nref.push_back(thisBoundaryCoordinates);
375}
376
377//_____________________________________________________________________________
378void RooMomentMorphND::Grid::addPdf(const RooAbsPdf &pdf, int bin_x, int bin_y)
379{
380 vector<int> thisBoundaries;
381 vector<double> thisBoundaryCoordinates;
382 thisBoundaries.push_back(bin_x);
383 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
384 thisBoundaries.push_back(bin_y);
385 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
386 _pdfList.add(pdf);
387 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
388 _nref.push_back(thisBoundaryCoordinates);
389}
390
391//_____________________________________________________________________________
392void RooMomentMorphND::Grid::addPdf(const RooAbsPdf &pdf, int bin_x, int bin_y, int bin_z)
393{
394 vector<int> thisBoundaries;
395 vector<double> thisBoundaryCoordinates;
396 thisBoundaries.push_back(bin_x);
397 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
398 thisBoundaries.push_back(bin_y);
399 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
400 thisBoundaries.push_back(bin_z);
401 thisBoundaryCoordinates.push_back(_grid[2]->array()[bin_z]);
402 _pdfList.add(pdf);
403 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
404 _nref.push_back(thisBoundaryCoordinates);
405}
406
407//_____________________________________________________________________________
408void RooMomentMorphND::Grid::addPdf(const RooAbsPdf &pdf, vector<int> bins)
409{
410 vector<double> thisBoundaryCoordinates;
411 int nBins = bins.size();
412 for (int i = 0; i < nBins; i++) {
413 thisBoundaryCoordinates.push_back(_grid[i]->array()[bins[i]]);
414 }
415 _pdfList.add(pdf);
416 _pdfMap[bins] = _pdfList.getSize() - 1;
417 _nref.push_back(thisBoundaryCoordinates);
418}
419
420//_____________________________________________________________________________
422{
423 CacheElem *cache = static_cast<CacheElem *>(_cacheMgr.getObj(0, (RooArgSet *)0));
424 if (cache) {
425 return cache;
426 }
427
428 int nObs = _obsList.getSize();
429 int nPdf = _referenceGrid._pdfList.getSize();
430
431 TIter pdfItr = _pdfList.createIterator();
432
433 RooAbsReal *null = 0;
434 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
435 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
436 vector<RooAbsReal *> myrms(nObs, null);
437 vector<RooAbsReal *> mypos(nObs, null);
438 vector<RooAbsReal *> slope(nPdf * nObs, null);
439 vector<RooAbsReal *> offsetrv(nPdf * nObs, null);
440 vector<RooAbsReal *> transVar(nPdf * nObs, null);
441 vector<RooAbsReal *> transPdf(nPdf, null);
442
443 RooArgSet ownedComps;
444 RooArgList fracl;
445
446 // fraction parameters
447 RooArgList coefList("coefList"); // fractions multiplied with input pdfs
448 RooArgList coefList2("coefList2"); // fractions multiplied with mean position of observable contribution
449 RooArgList coefList3("coefList3"); // fractions multiplied with rms position of observable contribution
450
451 for (int i = 0; i < 3 * nPdf; ++i) {
452 string fracName = Form("frac_%d", i);
453 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), /*value=*/1.); // to be set later
454
455 fracl.add(*frac);
456 if (i < nPdf)
457 coefList.add(*(RooRealVar *)(fracl.at(i)));
458 else if (i < 2 * nPdf)
459 coefList2.add(*(RooRealVar *)(fracl.at(i)));
460 else
461 coefList3.add(*(RooRealVar *)(fracl.at(i)));
462 ownedComps.add(*(RooRealVar *)(fracl.at(i)));
463 }
464
465 RooAddPdf *theSumPdf = 0;
466 string sumpdfName = Form("%s_sumpdf", GetName());
467
468 if (_useHorizMorph) {
469 // mean and sigma
470 RooArgList obsList(_obsList);
471 for (int i = 0; i < nPdf; ++i) {
472 for (int j = 0; j < nObs; ++j) {
473 RooAbsMoment *mom = nObs == 1 ? ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*obsList.at(j))
474 : ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*obsList.at(j), obsList);
475
478
479 sigmarv[sij(i, j)] = mom;
480 meanrv[sij(i, j)] = mom->mean();
481
482 ownedComps.add(*sigmarv[sij(i, j)]);
483 }
484 }
485
486 // slope and offset (to be set later, depend on nuisance parameters)
487 for (int j = 0; j < nObs; ++j) {
488 RooArgList meanList("meanList");
489 RooArgList rmsList("rmsList");
490 for (int i = 0; i < nPdf; ++i) {
491 meanList.add(*meanrv[sij(i, j)]);
492 rmsList.add(*sigmarv[sij(i, j)]);
493 }
494 string myrmsName = Form("%s_rms_%d", GetName(), j);
495 string myposName = Form("%s_pos_%d", GetName(), j);
496 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
497 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
498 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
499 }
500
501 // construction of unit pdfs
502 pdfItr.Reset();
503 RooAbsPdf *pdf;
504 RooArgList transPdfList;
505
506 for (int i = 0; i < nPdf; ++i) {
507 _obsItr->Reset();
508 RooRealVar *var;
509
510 pdf = (RooAbsPdf *)pdfItr.Next();
511 string pdfName = Form("pdf_%d", i);
512 RooCustomizer cust(*pdf, pdfName.c_str());
513
514 for (int j = 0; j < nObs; ++j) {
515 // slope and offset formulas
516 string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
517 string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
518
519 slope[sij(i, j)] =
520 new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[sij(i, j)], *myrms[j]));
521 offsetrv[sij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
522 RooArgList(*meanrv[sij(i, j)], *mypos[j], *slope[sij(i, j)]));
523 ownedComps.add(RooArgSet(*slope[sij(i, j)], *offsetrv[sij(i, j)]));
524
525 // linear transformations, so pdf can be renormalized easily
526 var = (RooRealVar *)(_obsItr->Next());
527 string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
528 transVar[sij(i, j)] = new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[sij(i, j)],
529 *offsetrv[sij(i, j)]);
530
531 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
532 // This will prevent the likelihood optimizers from erroneously declaring terms constant
533 transVar[sij(i, j)]->addServerList((RooAbsCollection &)_parList);
534
535 ownedComps.add(*transVar[sij(i, j)]);
536 cust.replaceArg(*var, *transVar[sij(i, j)]);
537 }
538 transPdf[i] = (RooAbsPdf *)cust.build();
539 transPdfList.add(*transPdf[i]);
540 ownedComps.add(*transPdf[i]);
541 }
542
543 // sum pdf
544 theSumPdf = new RooAddPdf(sumpdfName.c_str(), sumpdfName.c_str(), transPdfList, coefList);
545 } else {
546 theSumPdf = new RooAddPdf(sumpdfName.c_str(), sumpdfName.c_str(), _pdfList, coefList);
547 }
548
549 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
550 // This will prevent the likelihood optimizers from erroneously declaring terms constant
552 theSumPdf->addOwnedComponents(ownedComps);
553
554 // change tracker for fraction parameters
555 string trackerName = Form("%s_frac_tracker", GetName());
556 RooChangeTracker *tracker = new RooChangeTracker(trackerName.c_str(), trackerName.c_str(), _parList, kTRUE);
557
558 // Store it in the cache
559 cache = new CacheElem(*theSumPdf, *tracker, fracl);
560 _cacheMgr.setObj(0, 0, cache, 0);
561
562 cache->calculateFractions(*this, kFALSE);
563 return cache;
564}
565
566//_____________________________________________________________________________
568{
569 return RooArgList(*_sumPdf, *_tracker);
570}
571
572//_____________________________________________________________________________
574{
575 delete _sumPdf;
576 delete _tracker;
577}
578
579//_____________________________________________________________________________
581{
582 // Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
583 _curNormSet = set ? (RooArgSet *)set : (RooArgSet *)&_obsList;
584 return RooAbsPdf::getVal(set);
585}
586
587//_____________________________________________________________________________
589{
590 CacheElem *cache = getCache(nset ? nset : _curNormSet);
591
592 if (cache->_tracker->hasChanged(kTRUE)) {
593 cache->calculateFractions(*this, kFALSE); // verbose turned off
594 }
595 return cache->_sumPdf;
596}
597
598//_____________________________________________________________________________
600{
602
603 if (cache->_tracker->hasChanged(kTRUE)) {
604 cache->calculateFractions(*this, kFALSE); // verbose turned off
605 }
606
607 Double_t ret = cache->_sumPdf->getVal(_obsList.nset());
608
609 return ret;
610}
611
612//_____________________________________________________________________________
614{
615 return (RooRealVar *)(_frac.at(i));
616}
617
618//_____________________________________________________________________________
620{
621 return (RooRealVar *)(_frac.at(i));
622}
623
624//_____________________________________________________________________________
625// from http://stackoverflow.com/a/5097100/8747
626template <typename Iterator>
627inline bool next_combination(const Iterator first, Iterator k, const Iterator last)
628{
629 if ((first == last) || (first == k) || (last == k)) {
630 return false;
631 }
632 Iterator itr1 = first;
633 Iterator itr2 = last;
634 ++itr1;
635 if (last == itr1) {
636 return false;
637 }
638 itr1 = last;
639 --itr1;
640 itr1 = k;
641 --itr2;
642 while (first != itr1) {
643 if (*--itr1 < *itr2) {
644 Iterator j = k;
645 while (!(*itr1 < *j)) ++j;
646 iter_swap(itr1, j);
647 ++itr1;
648 ++j;
649 itr2 = k;
650 rotate(itr1, j, last);
651 while (last != j) {
652 ++j;
653 ++itr2;
654 }
655 rotate(k, itr2, last);
656 return true;
657 }
658 }
659 rotate(first, k, last);
660 return false;
661}
662
663//_____________________________________________________________________________
665{
666 int nPdf = self._pdfList.getSize();
667 int nPar = self._parList.getSize();
668
669 Double_t fracLinear(1.);
670 Double_t fracNonLinear(1.);
671
673 // Calculate the delta vector
674 vector<double> dm2;
675 for (int idim = 0; idim < nPar; idim++) {
676 Double_t delta = ((RooRealVar *)self._parList.at(idim))->getVal() - self._referenceGrid._nref[0][idim];
677 dm2.push_back(delta);
678 }
679
680 vector<vector<int>> powers;
681 for (int idim = 0; idim < nPar; idim++) {
682 vector<int> xtmp;
683 for (int ix = 0; ix < self._referenceGrid._nnuis[idim]; ix++) {
684 xtmp.push_back(ix);
685 }
686 powers.push_back(xtmp);
687 }
688
689 vector<vector<int>> output;
690 cartesian_product(output, powers);
691 int nCombs = output.size();
692
693 vector<double> deltavec(nPdf, 1.0);
694
695 int nperm = 0;
696 for (int i = 0; i < nCombs; i++) {
697 double tmpDm = 1.0;
698 for (int ix = 0; ix < nPar; ix++) {
699 Double_t delta = dm2[ix];
700 tmpDm *= TMath::Power(delta, static_cast<double>(output[i][ix]));
701 }
702 deltavec[nperm] = tmpDm;
703 nperm++;
704 }
705
706 double sumposfrac = 0.0;
707 for (int i = 0; i < nPdf; ++i) {
708 double ffrac = 0.0;
709
710 for (int j = 0; j < nPdf; ++j) {
711 ffrac += (*self._M)(j, i) * deltavec[j] * fracNonLinear;
712 }
713
714 if (ffrac >= 0) {
715 sumposfrac += ffrac;
716 }
717
718 // fractions for pdf
719 if (self._setting != NonLinearLinFractions) {
720 ((RooRealVar *)frac(i))->setVal(ffrac);
721 }
722
723 // fractions for rms and mean
724 ((RooRealVar *)frac(nPdf + i))->setVal(ffrac); // need to add up
725 ((RooRealVar *)frac(2 * nPdf + i))->setVal(ffrac); // need to add up
726
727 if (verbose) {
728 cout << "NonLinear fraction " << ffrac << endl;
729 frac(i)->Print();
730 frac(nPdf + i)->Print();
731 frac(2 * nPdf + i)->Print();
732 }
733 }
734
735 if (self._setting == NonLinearPosFractions) {
736 for (int i = 0; i < nPdf; ++i) {
737 if (((RooRealVar *)frac(i))->getVal() < 0)
738 ((RooRealVar *)frac(i))->setVal(0.);
739 ((RooRealVar *)frac(i))->setVal(((RooRealVar *)frac(i))->getVal() / sumposfrac);
740 }
741 }
742 }
743
744 if (self._setting == Linear || self._setting == NonLinearLinFractions) {
745 // loop over parList
746 self._parItr->Reset();
747
748 // zero all fractions
749 // for (int i = 0; i < 3*nPdf; ++i) {
750 for (int i = 0; i < nPdf; ++i) {
751 double initval = 0;
752 ((RooRealVar *)frac(i))->setVal(initval);
753 ((RooRealVar *)frac(nPdf + i))->setVal(initval);
754 ((RooRealVar *)frac(2 * nPdf + i))->setVal(initval);
755 }
756
757 vector<double> mtmp;
758
759 for (int j = 0; j < nPar; j++) {
760 RooRealVar *m = (RooRealVar *)(self._parItr->Next());
761 mtmp.push_back(m->getVal());
762 }
763
764 self.findShape(mtmp); // this sets _squareVec and _squareIdx quantities
765
766 int depth = TMath::Power(2, nPar);
767 vector<double> deltavec(depth, 1.0);
768
769 int nperm = 0;
770
771 vector<int> xtmp;
772 for (int ix = 0; ix < nPar; ix++) {
773 xtmp.push_back(ix);
774 }
775
776 for (int iperm = 1; iperm <= nPar; ++iperm) {
777 do {
778 double dtmp = mtmp[xtmp[0]] - self._squareVec[0][xtmp[0]];
779 for (int itmp = 1; itmp < iperm; ++itmp) {
780 dtmp *= mtmp[xtmp[itmp]] - self._squareVec[0][xtmp[itmp]];
781 }
782 deltavec[nperm + 1] = dtmp;
783 nperm++;
784 } while (next_combination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
785 }
786
787 Double_t origFrac1(0.), origFrac2(0.);
788 for (int i = 0; i < depth; ++i) {
789 double ffrac = 0.;
790 for (int j = 0; j < depth; ++j) {
791 ffrac += (*self._MSqr)(j, i) * deltavec[j] * fracLinear;
792 }
793
794 // set fractions for pdf
795 origFrac1 = ((RooRealVar *)frac(self._squareIdx[i]))->getVal(); // already set in case of smoothlinear
796 ((RooRealVar *)frac(self._squareIdx[i]))->setVal(origFrac1 + ffrac); // need to add up
797
798 // set fractions for rms and mean
799 if (self._setting != NonLinearLinFractions) {
800 origFrac2 =
801 ((RooRealVar *)frac(nPdf + self._squareIdx[i]))->getVal(); // already set in case of smoothlinear
802 ((RooRealVar *)frac(nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
803 ((RooRealVar *)frac(2 * nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
804 }
805
806 if (verbose) {
807 cout << "Linear fraction " << ffrac << endl;
808 frac(self._squareIdx[i])->Print();
809 frac(nPdf + self._squareIdx[i])->Print();
810 frac(2 * nPdf + self._squareIdx[i])->Print();
811 }
812 }
813 }
814}
815
816//_____________________________________________________________________________
817void RooMomentMorphND::findShape(const vector<double> &x) const
818{
819 int nPar = _parList.getSize();
820 int nRef = _referenceGrid._nref.size();
821
822 // Find hypercube enclosing the location to morph to
823 // bool isEnclosed = true;
824 // for (int i = 0; i < nPar; i++) {
825 // if (x[i] < _referenceGrid._grid[i]->lowBound())
826 // isEnclosed = false;
827 // if (x[i] > _referenceGrid._grid[i]->highBound())
828 // isEnclosed = false;
829 // }
830
831 // cout << "isEnclosed = " << isEnclosed << endl;
832
833 int depth = TMath::Power(2, nPar);
834
835 vector<vector<double>> boundaries(nPar);
836 for (int idim = 0; idim < nPar; idim++) {
837 int bin = _referenceGrid._grid[idim]->binNumber(x[idim]);
838 double lo = _referenceGrid._grid[idim]->binLow(bin);
839 double hi = _referenceGrid._grid[idim]->binHigh(bin);
840 boundaries[idim].push_back(lo);
841 boundaries[idim].push_back(hi);
842 }
843
844 vector<vector<double>> output;
845 cartesian_product(output, boundaries);
847
848 for (int isq = 0; isq < depth; isq++) {
849 for (int iref = 0; iref < nRef; iref++) {
850 if (_squareVec[isq] == _referenceGrid._nref[iref]) {
851 _squareIdx[isq] = iref;
852 break;
853 }
854 }
855 }
856
857 // cout << endl;
858
859 // for (int isq = 0; isq < _squareVec.size(); isq++) {
860 // cout << _squareIdx[isq];
861 // cout << " (";
862 // for (int isqq = 0; isqq < _squareVec[isq].size(); isqq++) {
863 // cout << _squareVec[isq][isqq] << ((isqq<_squareVec[isq].size()-1)?",":"");
864 // }
865 // cout << ") ";
866 // }
867
868 // construct transformation matrix for linear extrapolation
869 TMatrixD M(depth, depth);
870
871 vector<int> xtmp;
872 for (int ix = 0; ix < nPar; ix++) {
873 xtmp.push_back(ix);
874 }
875
876 for (int k = 0; k < depth; ++k) {
877 M(k, 0) = 1.0;
878
879 int nperm = 0;
880 vector<double> squareBase = _squareVec[0];
881
882 for (int iperm = 1; iperm <= nPar; ++iperm) {
883 do {
884 double dtmp = _squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
885 for (int itmp = 1; itmp < iperm; ++itmp) {
886 dtmp *= _squareVec[k][xtmp[itmp]] - squareBase[xtmp[itmp]];
887 }
888 M(k, nperm + 1) = dtmp;
889 nperm++;
890 } while (next_combination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
891 }
892 }
893
894 // M.Print();
895 (*_MSqr) = M.Invert();
896}
897
898//_____________________________________________________________________________
900{
901 if (allVars.getSize() == 1) {
902 RooAbsReal *temp = const_cast<RooMomentMorphND *>(this);
903 temp->specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator");
904 int nbins = ((RooRealVar *)allVars.first())->numBins();
905 temp->specialIntegratorConfig(kTRUE)->getConfigSection("RooBinIntegrator").setRealValue("numBins", nbins);
906 return true;
907 } else {
908 cout << "Currently BinIntegrator only knows how to deal with 1-d " << endl;
909 return false;
910 }
911 return false;
912}
#define d(i)
Definition RSha256.hxx:102
bool next_combination(const Iterator first, Iterator k, const Iterator last)
void cartesian_product(vector< vector< T > > &out, vector< vector< T > > &in)
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
#define hi
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
friend class RooArgSet
Definition RooAbsArg.h:642
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
void addServerList(RooAbsCollection &serverList, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE)
Register a list of RooAbsArg as servers to us by calling addServer() for each arg in the list.
void setLocalNoDirtyInhibit(Bool_t flag) const
Definition RooAbsArg.h:720
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
RooAbsArg * first() const
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
RooAbsReal * mean()
const RooArgSet * nset() const
Definition RooAbsProxy.h:45
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
friend class RooAddPdf
Definition RooAbsReal.h:558
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooAddition 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:35
Class RooBinning is an implements RooAbsBinning in terms of an array of boundary values,...
Definition RooBinning.h:28
virtual Int_t numBoundaries() const
Definition RooBinning.h:38
virtual Double_t * array() const
Return array of boundary values.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
virtual Bool_t setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
Bool_t hasChanged(Bool_t clearState)
Returns true if state has changed since last call with clearState=kTRUE.
RooConstVar represent a constant real-valued object.
Definition RooConstVar.h:26
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
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...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
virtual RooArgList containedArgs(Action)
void calculateFractions(const RooMomentMorphND &self, Bool_t verbose=kTRUE) const
void addBinning(const RooAbsBinning &binning)
std::vector< std::vector< double > > _nref
std::vector< RooAbsBinning * > _grid
void addPdf(const RooAbsPdf &pdf, int bin_x)
std::vector< int > _nnuis
virtual Double_t getVal(const RooArgSet *set=0) const
std::vector< int > _squareIdx
Grid _referenceGrid
Do not persist.
CacheElem * getCache(const RooArgSet *nset) const
RooObjCacheManager _cacheMgr
! Transient cache manager
RooArgSet * _curNormSet
Bool_t setBinIntegrator(RooArgSet &allVars)
std::vector< std::vector< double > > _squareVec
void initializeObservables(const RooArgList &obsList)
int sij(const int &i, const int &j) const
void findShape(const std::vector< double > &x) const
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
TIterator * _obsItr
Do not persist.
RooAbsPdf * sumPdf(const RooArgSet *nset)
void initializeParameters(const RooArgList &parList)
RooListProxy _parList
RooListProxy _pdfList
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooCategory & method1D()
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
TObject * Next()
void Reset()
Iterator abstract base class.
Definition TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Int_t GetNrows() const
Definition TVectorT.h:75
Element * GetMatrixArray()
Definition TVectorT.h:78
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:685
Definition first.py:1
vector< T >::const_iterator begin
vector< T >::const_iterator end
vector< T >::const_iterator me
auto * m
Definition textangle.C:8
static void output(int code)
Definition gifencode.c:226