Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorphFuncND.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
10#include "RooAbsCategory.h"
11#include "RooAddPdf.h"
12#include "RooRealIntegral.h"
13#include "RooRealConstant.h"
14#include "RooRealVar.h"
15#include "RooFormulaVar.h"
16#include "RooCustomizer.h"
17#include "RooRealSumFunc.h"
18#include "RooAddition.h"
19#include "RooAbsMoment.h"
20#include "RooMoment.h"
21#include "RooLinearVar.h"
22#include "RooChangeTracker.h"
23#include "RooNumIntConfig.h"
24#include "RooHistPdf.h"
25
27
28#include "TMath.h"
29#include "TVector.h"
30#include "TMap.h"
31
32#include <map>
33#include <algorithm>
34
35using namespace std;
36
39
40//_____________________________________________________________________________
41RooMomentMorphFuncND::RooMomentMorphFuncND() : _cacheMgr(this, 10, true, true), _setting(RooMomentMorphFuncND::Linear), _useHorizMorph(true)
42{
44}
45
46//_____________________________________________________________________________
47RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, const RooArgList &parList, const RooArgList &obsList,
48 const Grid2 &referenceGrid, const Setting &setting)
50 _cacheMgr(this, 10, true, true),
51 _parList("parList", "List of morph parameters", this),
52 _obsList("obsList", "List of observables", this),
53 _referenceGrid(referenceGrid),
54 _pdfList("pdfList", "List of pdfs", this),
55 _setting(setting),
56 _useHorizMorph(true)
57{
58 // morph parameters
59 initializeParameters(parList);
60
61 // observables
62 initializeObservables(obsList);
63
65
66 // general initialization
67 initialize();
68
70}
71
72//_____________________________________________________________________________
73RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
74 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
76 _cacheMgr(this, 10, true, true),
77 _parList("parList", "List of morph parameters", this),
78 _obsList("obsList", "List of observables", this),
79 _pdfList("pdfList", "List of pdfs", this),
80 _setting(setting),
81 _useHorizMorph(true)
82{
83 // make reference grid
84 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
86
87 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
88 for (int j = 0; j < grid.numBoundaries(); ++j) {
89 if (mrefpoints[i] == grid.array()[j]) {
90 _referenceGrid.addPdf(*(Base_t *)pdfList.at(i), j);
91 break;
92 }
93 }
94 }
95
97
98 // morph parameters
99 RooArgList parList;
100 parList.add(_m);
101 initializeParameters(parList);
102
103 // observables
104 initializeObservables(varList);
105
106 // general initialization
107 initialize();
108
110}
111
112//_____________________________________________________________________________
113RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
114 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
116 _cacheMgr(this, 10, true, true),
117 _parList("parList", "List of morph parameters", this),
118 _obsList("obsList", "List of observables", this),
119 _pdfList("pdfList", "List of pdfs", this),
120 _setting(setting),
121 _useHorizMorph(true)
122{
123 // make reference grid
124 TVectorD mrefpoints(mrefList.getSize());
125 Int_t i = 0;
126 for (auto *mref : mrefList) {
127 if (!dynamic_cast<RooAbsReal *>(mref)) {
128 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
129 << " is not of type RooAbsReal" << endl;
130 throw string("RooMomentMorphFuncND::ctor() ERROR mref is not of type RooAbsReal");
131 }
132 if (!dynamic_cast<RooConstVar *>(mref)) {
133 coutW(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") WARNING mref point " << i
134 << " is not a constant, taking a snapshot of its value" << endl;
135 }
136 mrefpoints[i] = static_cast<RooAbsReal *>(mref)->getVal();
137 i++;
138 }
139
140 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
142 for (i = 0; i < mrefpoints.GetNrows(); ++i) {
143 for (int j = 0; j < grid.numBoundaries(); ++j) {
144 if (mrefpoints[i] == grid.array()[j]) {
145 _referenceGrid.addPdf(static_cast<Base_t &>(pdfList[i]), j);
146 break;
147 }
148 }
149 }
150
152
153 // morph parameters
154 RooArgList parList;
155 parList.add(_m);
156 initializeParameters(parList);
157
158 // observables
159 initializeObservables(varList);
160
161 // general initialization
162 initialize();
163
165}
166
167//_____________________________________________________________________________
170 _cacheMgr(other._cacheMgr, this),
171 _parList("parList", this, other._parList),
172 _obsList("obsList", this, other._obsList),
173 _referenceGrid(other._referenceGrid),
174 _pdfList("pdfList", this, other._pdfList),
175 _setting(other._setting),
176 _useHorizMorph(other._useHorizMorph)
177{
178 // general initialization
179 initialize();
180
182}
183
184//_____________________________________________________________________________
186{
188}
189
190//_____________________________________________________________________________
192{
193 for (auto *par : parList) {
194 if (!dynamic_cast<RooAbsReal *>(par)) {
195 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: parameter " << par->GetName()
196 << " is not of type RooAbsReal" << endl;
197 throw string("RooMomentMorphFuncND::initializeParameters() ERROR parameter is not of type RooAbsReal");
198 }
199 _parList.add(*par);
200 }
201}
202
203//_____________________________________________________________________________
205{
206 for (auto *var : obsList) {
207 if (!dynamic_cast<RooAbsReal *>(var)) {
208 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
209 << " is not of type RooAbsReal" << endl;
210 throw string("RooMomentMorphFuncND::initializeObservables() ERROR variable is not of type RooAbsReal");
211 }
212 _obsList.add(*var);
213 }
214}
215
216//_____________________________________________________________________________
218{
219 for (vector<RooAbsBinning *>::iterator itr = _referenceGrid._grid.begin(); itr != _referenceGrid._grid.end();
220 ++itr) {
221 _referenceGrid._nnuis.push_back((*itr)->numBins() + 1);
222 }
223
224 int nPar = _parList.getSize();
225 int nDim = _referenceGrid._grid.size();
226 int nPdf = _referenceGrid._pdfList.getSize();
227 int nRef = _referenceGrid._nref.size();
228 int depth = TMath::Power(2, nPar);
229
230 if (nPar != nDim) {
231 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPar != nDim"
232 << ": " << nPar << " !=" << nDim << endl;
233 assert(0);
234 }
235
236 if (nPdf != nRef) {
237 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPdf != nRef"
238 << ": " << nPdf << " !=" << nRef << endl;
239 assert(0);
240 }
241
242 // Transformation matrix for NonLinear settings
243 _M = std::make_unique<TMatrixD>(nPdf, nPdf);
244 _MSqr = std::make_unique<TMatrixD>(depth, depth);
246 TMatrixD M(nPdf, nPdf);
247
248 vector<vector<double>> dm(nPdf);
249 for (int k = 0; k < nPdf; ++k) {
250 vector<double> dm2;
251 for (int idim = 0; idim < nPar; idim++) {
252 double delta = _referenceGrid._nref[k][idim] - _referenceGrid._nref[0][idim];
253 dm2.push_back(delta);
254 }
255 dm[k] = dm2;
256 }
257
258 vector<vector<int>> powers;
259 for (int idim = 0; idim < nPar; idim++) {
260 vector<int> xtmp;
261 xtmp.reserve(_referenceGrid._nnuis[idim]);
262 for (int ix = 0; ix < _referenceGrid._nnuis[idim]; ix++) {
263 xtmp.push_back(ix);
264 }
265 powers.push_back(xtmp);
266 }
267
268 vector<vector<int>> output;
270 int nCombs = output.size();
271
272 for (int k = 0; k < nPdf; ++k) {
273 int nperm = 0;
274 for (int i = 0; i < nCombs; i++) {
275 double tmpDm = 1.0;
276 for (int ix = 0; ix < nPar; ix++) {
277 double delta = dm[k][ix];
278 tmpDm *= TMath::Power(delta, static_cast<double>(output[i][ix]));
279 }
280 M(k, nperm) = tmpDm;
281 nperm++;
282 }
283 }
284
285 // M.Print();
286 (*_M) = M.Invert();
287 }
288
289 // Resize transformation vectors
290 _squareVec.resize(TMath::Power(2, nPar));
291 _squareIdx.resize(TMath::Power(2, nPar));
292}
293
294//_____________________________________________________________________________
296 : _pdfList(other._pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
297{
298 for (unsigned int i = 0; i < other._grid.size(); i++) {
299 _grid.push_back(other._grid[i]->clone());
300 }
301}
302
303//_____________________________________________________________________________
305{
306 for (RooAbsBinning *binning : _grid) {
307 delete binning;
308 }
309}
310
311//_____________________________________________________________________________
313{
314 vector<int> thisBoundaries;
315 vector<double> thisBoundaryCoordinates;
316 thisBoundaries.push_back(bin_x);
317 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
318 _pdfList.add(pdf);
319 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
320 _nref.push_back(thisBoundaryCoordinates);
321}
322
323//_____________________________________________________________________________
325{
326 vector<int> thisBoundaries;
327 vector<double> thisBoundaryCoordinates;
328 thisBoundaries.push_back(bin_x);
329 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
330 thisBoundaries.push_back(bin_y);
331 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
332 _pdfList.add(pdf);
333 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
334 _nref.push_back(thisBoundaryCoordinates);
335}
336
337//_____________________________________________________________________________
338void RooMomentMorphFuncND::Grid2::addPdf(const RooMomentMorphFuncND::Base_t &pdf, int bin_x, int bin_y, int bin_z)
339{
340 vector<int> thisBoundaries;
341 vector<double> thisBoundaryCoordinates;
342 thisBoundaries.push_back(bin_x);
343 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
344 thisBoundaries.push_back(bin_y);
345 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
346 thisBoundaries.push_back(bin_z);
347 thisBoundaryCoordinates.push_back(_grid[2]->array()[bin_z]);
348 _pdfList.add(pdf);
349 _pdfMap[thisBoundaries] = _pdfList.getSize() - 1;
350 _nref.push_back(thisBoundaryCoordinates);
351}
352
353//_____________________________________________________________________________
355{
356 vector<double> thisBoundaryCoordinates;
357 int nBins = bins.size();
358 thisBoundaryCoordinates.reserve(nBins);
359 for (int i = 0; i < nBins; i++) {
360 thisBoundaryCoordinates.push_back(_grid[i]->array()[bins[i]]);
361 }
362 _pdfList.add(pdf);
363 _pdfMap[bins] = _pdfList.getSize() - 1;
364 _nref.push_back(thisBoundaryCoordinates);
365}
366
367//_____________________________________________________________________________
369{
370 auto cache = static_cast<CacheElem *>(_cacheMgr.getObj(nullptr, static_cast<RooArgSet const*>(nullptr)));
371 if (cache) {
372 return cache;
373 }
374
375 int nObs = _obsList.getSize();
376 int nPdf = _referenceGrid._pdfList.getSize();
377
378 RooAbsReal *null = nullptr;
379 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
380 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
381 vector<RooAbsReal *> myrms(nObs, null);
382 vector<RooAbsReal *> mypos(nObs, null);
383 vector<RooAbsReal *> slope(nPdf * nObs, null);
384 vector<RooAbsReal *> offsets(nPdf * nObs, null);
385 vector<RooAbsReal *> transVar(nPdf * nObs, null);
386 vector<RooAbsReal *> transPdf(nPdf, null);
387
388 RooArgSet ownedComps;
389 RooArgList fracl;
390
391 // fraction parameters
392 RooArgList coefList("coefList"); // fractions multiplied with input pdfs
393 RooArgList coefList2("coefList2"); // fractions multiplied with mean position of observable contribution
394 RooArgList coefList3("coefList3"); // fractions multiplied with rms position of observable contribution
395
396 for (int i = 0; i < 3 * nPdf; ++i) {
397 string fracName = Form("frac_%d", i);
398 double initval = _isPdfMode ? 1.0 : 0.0;
399 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), initval); // to be set later
400
401 fracl.add(*frac);
402 if (i < nPdf)
403 coefList.add(*(RooRealVar *)(fracl.at(i)));
404 else if (i < 2 * nPdf)
405 coefList2.add(*(RooRealVar *)(fracl.at(i)));
406 else
407 coefList3.add(*(RooRealVar *)(fracl.at(i)));
408 ownedComps.add(*(RooRealVar *)(fracl.at(i)));
409 }
410
411 RooAbsReal *theSum = nullptr;
412 string sumName = Form("%s_sum", GetName());
413
414 RooArgList transPdfList;
415 if (_useHorizMorph) {
416 // mean and sigma
417 RooArgList obsList(_obsList);
418 for (int i = 0; i < nPdf; ++i) {
419 for (int j = 0; j < nObs; ++j) {
420 RooAbsMoment *mom = nObs == 1 ? ((Base_t *)_pdfList.at(i))->sigma((RooRealVar &)*obsList.at(j))
421 : ((Base_t *)_pdfList.at(i))->sigma((RooRealVar &)*obsList.at(j), obsList);
422
423 mom->setLocalNoDirtyInhibit(true);
424 mom->mean()->setLocalNoDirtyInhibit(true);
425
426 sigmarv[sij(i, j)] = mom;
427 meanrv[sij(i, j)] = mom->mean();
428
429 ownedComps.add(*sigmarv[sij(i, j)]);
430 }
431 }
432
433 // slope and offset (to be set later, depend on nuisance parameters)
434 for (int j = 0; j < nObs; ++j) {
435 RooArgList meanList("meanList");
436 RooArgList rmsList("rmsList");
437 for (int i = 0; i < nPdf; ++i) {
438 meanList.add(*meanrv[sij(i, j)]);
439 rmsList.add(*sigmarv[sij(i, j)]);
440 }
441 string myrmsName = Form("%s_rms_%d", GetName(), j);
442 string myposName = Form("%s_pos_%d", GetName(), j);
443 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
444 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
445 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
446 }
447
448 // construction of unit pdfs
449
450 Int_t i = 0;
451 for (auto const *pdf : static_range_cast<Base_t *>(_pdfList)) {
452
453 string pdfName = Form("pdf_%d", i);
454 RooCustomizer cust(*pdf, pdfName.c_str());
455
456 Int_t j = 0;
457 for (auto *var : static_range_cast<RooRealVar *>(obsList)) {
458 // slope and offset formulas
459 string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
460 string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
461
462 slope[sij(i, j)] =
463 new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[sij(i, j)], *myrms[j]));
464 offsets[sij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
465 RooArgList(*meanrv[sij(i, j)], *mypos[j], *slope[sij(i, j)]));
466 ownedComps.add(RooArgSet(*slope[sij(i, j)], *offsets[sij(i, j)]));
467
468 // linear transformations, so pdf can be renormalized easily
469 string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
470 transVar[sij(i, j)] = new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[sij(i, j)],
471 *offsets[sij(i, j)]);
472
473 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
474 // This will prevent the likelihood optimizers from erroneously declaring terms constant
475 transVar[sij(i, j)]->addServerList((RooAbsCollection &)_parList);
476
477 ownedComps.add(*transVar[sij(i, j)]);
478 cust.replaceArg(*var, *transVar[sij(i, j)]);
479 ++j;
480 }
481 transPdf[i] = static_cast<Base_t *>(cust.build());
482 transPdfList.add(*transPdf[i]);
483 ownedComps.add(*transPdf[i]);
484 ++i;
485 }
486 }
487
488 // sum pdf
489 RooArgList const &pdfList = _useHorizMorph ? transPdfList : static_cast<RooArgList const &>(_pdfList);
490 if (_isPdfMode) {
491 theSum = new RooAddPdf(sumName.c_str(), sumName.c_str(), pdfList, coefList);
492 } else {
493 theSum = new RooRealSumFunc(sumName.c_str(), sumName.c_str(), pdfList, coefList);
494 }
495
496 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
497 // This will prevent the likelihood optimizers from erroneously declaring terms constant
499 theSum->addOwnedComponents(ownedComps);
500
501 // change tracker for fraction parameters
502 string trackerName = Form("%s_frac_tracker", GetName());
503 RooChangeTracker *tracker = new RooChangeTracker(trackerName.c_str(), trackerName.c_str(), _parList, true);
504
505 // Store it in the cache
506 cache = new CacheElem(*theSum, *tracker, fracl);
507 _cacheMgr.setObj(nullptr, nullptr, cache, nullptr);
508
509 return cache;
510}
511
512//_____________________________________________________________________________
514{
515 return RooArgList(*_sum, *_tracker);
516}
517
518//_____________________________________________________________________________
520{
521 delete _sum;
522 delete _tracker;
523}
524
525//_____________________________________________________________________________
527{
528 // Special version of getValV() overrides Base_t::getValV() to save value of current normalization set
529 _curNormSet = set ? (RooArgSet *)set : (RooArgSet *)&_obsList;
530 return Base_t::getValV(set);
531}
532
533//_____________________________________________________________________________
535{
536 CacheElem *cache = getCache(nset ? nset : _curNormSet);
537
538 if (cache->_tracker->hasChanged(true)) {
539 cache->calculateFractions(*this, false); // verbose turned off
540 }
541 return cache->_sum;
542}
543
544//_____________________________________________________________________________
546{
548
549 if (cache->_tracker->hasChanged(true)) {
550 cache->calculateFractions(*this, false); // verbose turned off
551 }
552
553 double ret = cache->_sum->getVal(_obsList.nset());
554
555 return ret;
556}
557
558//_____________________________________________________________________________
560{
561 return (RooRealVar *)(_frac.at(i));
562}
563
564//_____________________________________________________________________________
566{
567 return (RooRealVar *)(_frac.at(i));
568}
569
570//_____________________________________________________________________________
572{
573 int nPdf = self._pdfList.getSize();
574 int nPar = self._parList.getSize();
575
576 double fracLinear(1.);
577 double fracNonLinear(1.);
578
580 // Calculate the delta vector
581 vector<double> dm2;
582 for (int idim = 0; idim < nPar; idim++) {
583 double delta = ((RooRealVar *)self._parList.at(idim))->getVal() - self._referenceGrid._nref[0][idim];
584 dm2.push_back(delta);
585 }
586
587 vector<vector<int>> powers;
588 for (int idim = 0; idim < nPar; idim++) {
589 vector<int> xtmp;
590 xtmp.reserve(self._referenceGrid._nnuis[idim]);
591 for (int ix = 0; ix < self._referenceGrid._nnuis[idim]; ix++) {
592 xtmp.push_back(ix);
593 }
594 powers.push_back(xtmp);
595 }
596
597 vector<vector<int>> output;
599 int nCombs = output.size();
600
601 vector<double> deltavec(nPdf, 1.0);
602
603 int nperm = 0;
604 for (int i = 0; i < nCombs; i++) {
605 double tmpDm = 1.0;
606 for (int ix = 0; ix < nPar; ix++) {
607 double delta = dm2[ix];
608 tmpDm *= TMath::Power(delta, static_cast<double>(output[i][ix]));
609 }
610 deltavec[nperm] = tmpDm;
611 nperm++;
612 }
613
614 double sumposfrac = 0.0;
615 for (int i = 0; i < nPdf; ++i) {
616 double ffrac = 0.0;
617
618 for (int j = 0; j < nPdf; ++j) {
619 ffrac += (*self._M)(j, i) * deltavec[j] * fracNonLinear;
620 }
621
622 if (ffrac >= 0) {
623 sumposfrac += ffrac;
624 }
625
626 // fractions for pdf
627 if (self._setting != NonLinearLinFractions) {
628 ((RooRealVar *)frac(i))->setVal(ffrac);
629 }
630
631 // fractions for rms and mean
632 ((RooRealVar *)frac(nPdf + i))->setVal(ffrac); // need to add up
633 ((RooRealVar *)frac(2 * nPdf + i))->setVal(ffrac); // need to add up
634
635 if (verbose) {
636 cout << "NonLinear fraction " << ffrac << endl;
637 frac(i)->Print();
638 frac(nPdf + i)->Print();
639 frac(2 * nPdf + i)->Print();
640 }
641 }
642
643 if (self._setting == NonLinearPosFractions) {
644 for (int i = 0; i < nPdf; ++i) {
645 if (((RooRealVar *)frac(i))->getVal() < 0)
646 ((RooRealVar *)frac(i))->setVal(0.);
647 ((RooRealVar *)frac(i))->setVal(((RooRealVar *)frac(i))->getVal() / sumposfrac);
648 }
649 }
650 }
651
652 if (self._setting == Linear || self._setting == NonLinearLinFractions) {
653 // zero all fractions
654 // for (int i = 0; i < 3*nPdf; ++i) {
655 for (int i = 0; i < nPdf; ++i) {
656 double initval = 0;
657 ((RooRealVar *)frac(i))->setVal(initval);
658 ((RooRealVar *)frac(nPdf + i))->setVal(initval);
659 ((RooRealVar *)frac(2 * nPdf + i))->setVal(initval);
660 }
661
662 std::vector<double> mtmp;
663
664 // loop over parList
665 for (auto *m : static_range_cast<RooRealVar *>(self._parList)) {
666 mtmp.push_back(m->getVal());
667 }
668
669 self.findShape(mtmp); // this sets _squareVec and _squareIdx quantities
670
671 int depth = TMath::Power(2, nPar);
672 vector<double> deltavec(depth, 1.0);
673
674 int nperm = 0;
675
676 vector<int> xtmp;
677 xtmp.reserve(nPar);
678 for (int ix = 0; ix < nPar; ix++) {
679 xtmp.push_back(ix);
680 }
681
682 for (int iperm = 1; iperm <= nPar; ++iperm) {
683 do {
684 double dtmp = mtmp[xtmp[0]] - self._squareVec[0][xtmp[0]];
685 for (int itmp = 1; itmp < iperm; ++itmp) {
686 dtmp *= mtmp[xtmp[itmp]] - self._squareVec[0][xtmp[itmp]];
687 }
688 deltavec[nperm + 1] = dtmp;
689 nperm++;
690 } while (RooFit::Detail::nextCombination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
691 }
692
693 double origFrac1(0.), origFrac2(0.);
694 for (int i = 0; i < depth; ++i) {
695 double ffrac = 0.;
696 for (int j = 0; j < depth; ++j) {
697 ffrac += (*self._MSqr)(j, i) * deltavec[j] * fracLinear;
698 }
699
700 // set fractions for pdf
701 origFrac1 = ((RooRealVar *)frac(self._squareIdx[i]))->getVal(); // already set in case of smoothlinear
702 ((RooRealVar *)frac(self._squareIdx[i]))->setVal(origFrac1 + ffrac); // need to add up
703
704 // set fractions for rms and mean
705 if (self._setting != NonLinearLinFractions) {
706 origFrac2 =
707 ((RooRealVar *)frac(nPdf + self._squareIdx[i]))->getVal(); // already set in case of smoothlinear
708 ((RooRealVar *)frac(nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
709 ((RooRealVar *)frac(2 * nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
710 }
711
712 if (verbose) {
713 cout << "Linear fraction " << ffrac << endl;
714 frac(self._squareIdx[i])->Print();
715 frac(nPdf + self._squareIdx[i])->Print();
716 frac(2 * nPdf + self._squareIdx[i])->Print();
717 }
718 }
719 }
720}
721
722//_____________________________________________________________________________
723void RooMomentMorphFuncND::findShape(const vector<double> &x) const
724{
725 int nPar = _parList.getSize();
726 int nRef = _referenceGrid._nref.size();
727
728 // Find hypercube enclosing the location to morph to
729 // bool isEnclosed = true;
730 // for (int i = 0; i < nPar; i++) {
731 // if (x[i] < _referenceGrid._grid[i]->lowBound())
732 // isEnclosed = false;
733 // if (x[i] > _referenceGrid._grid[i]->highBound())
734 // isEnclosed = false;
735 // }
736
737 // cout << "isEnclosed = " << isEnclosed << endl;
738
739 int depth = TMath::Power(2, nPar);
740
741 vector<vector<double>> boundaries(nPar);
742 for (int idim = 0; idim < nPar; idim++) {
743 int bin = _referenceGrid._grid[idim]->binNumber(x[idim]);
744 double lo = _referenceGrid._grid[idim]->binLow(bin);
745 double hi = _referenceGrid._grid[idim]->binHigh(bin);
746 boundaries[idim].push_back(lo);
747 boundaries[idim].push_back(hi);
748 }
749
750 vector<vector<double>> output;
753
754 for (int isq = 0; isq < depth; isq++) {
755 for (int iref = 0; iref < nRef; iref++) {
756 if (_squareVec[isq] == _referenceGrid._nref[iref]) {
757 _squareIdx[isq] = iref;
758 break;
759 }
760 }
761 }
762
763 // cout << endl;
764
765 // for (int isq = 0; isq < _squareVec.size(); isq++) {
766 // cout << _squareIdx[isq];
767 // cout << " (";
768 // for (int isqq = 0; isqq < _squareVec[isq].size(); isqq++) {
769 // cout << _squareVec[isq][isqq] << ((isqq<_squareVec[isq].size()-1)?",":"");
770 // }
771 // cout << ") ";
772 // }
773
774 // construct transformation matrix for linear extrapolation
775 TMatrixD M(depth, depth);
776
777 vector<int> xtmp;
778 xtmp.reserve(nPar);
779 for (int ix = 0; ix < nPar; ix++) {
780 xtmp.push_back(ix);
781 }
782
783 for (int k = 0; k < depth; ++k) {
784 M(k, 0) = 1.0;
785
786 int nperm = 0;
787 vector<double> squareBase = _squareVec[0];
788
789 for (int iperm = 1; iperm <= nPar; ++iperm) {
790 do {
791 double dtmp = _squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
792 for (int itmp = 1; itmp < iperm; ++itmp) {
793 dtmp *= _squareVec[k][xtmp[itmp]] - squareBase[xtmp[itmp]];
794 }
795 M(k, nperm + 1) = dtmp;
796 nperm++;
797 } while (RooFit::Detail::nextCombination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
798 }
799 }
800
801 // M.Print();
802 (*_MSqr) = M.Invert();
803}
804
805//_____________________________________________________________________________
807{
808 if (allVars.getSize() == 1) {
809 RooAbsReal *temp = const_cast<RooMomentMorphFuncND *>(this);
810 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator");
811 int nbins = ((RooRealVar *)allVars.first())->numBins();
812 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins", nbins);
813 return true;
814 } else {
815 cout << "Currently BinIntegrator only knows how to deal with 1-d " << endl;
816 return false;
817 }
818 return false;
819}
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
#define hi
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
void setLocalNoDirtyInhibit(bool flag) const
Definition RooAbsArg.h:702
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
void addServerList(RooAbsCollection &serverList, bool valueProp=true, bool shapeProp=false)
Register a list of RooAbsArg as servers to us by calling addServer() for each arg in the list.
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
Abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * first() const
bool setRealValue(const char *name, double newVal=0.0, bool verbose=false)
Set value of a RooAbsRealLValue stored in set with given name to newVal No error messages are printed...
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
RooAbsReal * mean()
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:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual double getValV(const RooArgSet *normalisationSet=nullptr) const
Return value of object.
friend class RooRealSumFunc
Definition RooAbsReal.h:410
friend class RooAddPdf
Definition RooAbsReal.h:412
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
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:55
Class RooBinning is an implements RooAbsBinning in terms of an array of boundary values,...
Definition RooBinning.h:27
double * array() const override
Return array of boundary values.
Int_t numBoundaries() const override
Return the number boundaries.
Definition RooBinning.h:38
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.
Meta object that tracks value changes in a given set of RooAbsArgs by registering itself as value cli...
bool hasChanged(bool clearState)
Returns true if state has changed since last call with clearState=true.
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...
RooConstVar represent 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 replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurrence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, bool verbose=false)
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...
void calculateFractions(const RooMomentMorphFuncND &self, bool verbose=true) const
RooArgList containedArgs(Action) override
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.
void initializeParameters(const RooArgList &parList)
RooAbsReal * sumFunc(const RooArgSet *nset)
CacheElem * getCache(const RooArgSet *nset) const
std::unique_ptr< TMatrixD > _M
RooArgSet * _curNormSet
! Transient cache manager
void initializeObservables(const RooArgList &obsList)
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::vector< int > _squareIdx
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:37
TMatrixT< Element > & Invert(Double_t *det=nullptr)
Invert the matrix and calculate its determinant.
const char * GetName() const override
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
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
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
TMarker m
Definition textangle.C:8
static void output()