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