Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooMomentMorphFuncND.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2023, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11/** \class RooMomentMorphFuncND
12 \ingroup Roofit
13
14**/
15
17
18#include <RooAbsCategory.h>
19#include <RooAbsMoment.h>
20#include <RooAddPdf.h>
21#include <RooAddition.h>
22#include <RooChangeTracker.h>
23#include <RooConstVar.h>
24#include <RooCustomizer.h>
25#include <RooFormulaVar.h>
26#include <RooLinearVar.h>
27#include <RooMoment.h>
28#include <RooNumIntConfig.h>
29#include <RooRealSumFunc.h>
30#include <RooRealVar.h>
31
33
34#include <Riostream.h>
35
36#include <TMap.h>
37#include <TMath.h>
38#include <TVector.h>
39
40#include <algorithm>
41#include <map>
42
43using std::string, std::vector;
44
45
46//_____________________________________________________________________________
47RooMomentMorphFuncND::RooMomentMorphFuncND() : _cacheMgr(this, 10, true, true), _setting(RooMomentMorphFuncND::Linear), _useHorizMorph(true)
48{
50}
51
52//_____________________________________________________________________________
56 _cacheMgr(this, 10, true, true),
57 _parList("parList", "List of morph parameters", this),
58 _obsList("obsList", "List of observables", this),
59 _referenceGrid(referenceGrid),
60 _pdfList("pdfList", "List of pdfs", this),
61 _setting(setting),
62 _useHorizMorph(true)
63{
64 // morph parameters
66
67 // observables
69
71
72 // general initialization
73 initialize();
74
76}
77
78//_____________________________________________________________________________
80 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
82 _cacheMgr(this, 10, true, true),
83 _parList("parList", "List of morph parameters", this),
84 _obsList("obsList", "List of observables", this),
85 _pdfList("pdfList", "List of pdfs", this),
86 _setting(setting),
87 _useHorizMorph(true)
88{
89 // make reference grid
90 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
92
93 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
94 for (int j = 0; j < grid.numBoundaries(); ++j) {
95 if (mrefpoints[i] == grid.array()[j]) {
96 _referenceGrid.addPdf(*static_cast<Base_t *>(pdfList.at(i)), j);
97 break;
98 }
99 }
100 }
101
103
104 // morph parameters
106 parList.add(_m);
108
109 // observables
111
112 // general initialization
113 initialize();
114
116}
117
118//_____________________________________________________________________________
120 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
122 _cacheMgr(this, 10, true, true),
123 _parList("parList", "List of morph parameters", this),
124 _obsList("obsList", "List of observables", this),
125 _pdfList("pdfList", "List of pdfs", this),
126 _setting(setting),
127 _useHorizMorph(true)
128{
129 // make reference grid
131 Int_t i = 0;
132 for (auto *mref : mrefList) {
133 if (!dynamic_cast<RooAbsReal *>(mref)) {
134 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
135 << " is not of type RooAbsReal" << std::endl;
136 throw string("RooMomentMorphFuncND::ctor() ERROR mref is not of type RooAbsReal");
137 }
138 if (!dynamic_cast<RooConstVar *>(mref)) {
139 coutW(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") WARNING mref point " << i
140 << " is not a constant, taking a snapshot of its value" << std::endl;
141 }
142 mrefpoints[i] = static_cast<RooAbsReal *>(mref)->getVal();
143 i++;
144 }
145
146 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
148 for (i = 0; i < mrefpoints.GetNrows(); ++i) {
149 for (int j = 0; j < grid.numBoundaries(); ++j) {
150 if (mrefpoints[i] == grid.array()[j]) {
151 _referenceGrid.addPdf(static_cast<Base_t &>(pdfList[i]), j);
152 break;
153 }
154 }
155 }
156
158
159 // morph parameters
161 parList.add(_m);
163
164 // observables
166
167 // general initialization
168 initialize();
169
171}
172
173//_____________________________________________________________________________
176 _cacheMgr(other._cacheMgr, this),
177 _parList("parList", this, other._parList),
178 _obsList("obsList", this, other._obsList),
179 _referenceGrid(other._referenceGrid),
180 _pdfList("pdfList", this, other._pdfList),
181 _setting(other._setting),
182 _useHorizMorph(other._useHorizMorph),
183 _isPdfMode{other._isPdfMode}
184{
185 // general initialization
186 initialize();
187
189}
190
191//_____________________________________________________________________________
196
197//_____________________________________________________________________________
199{
200 for (vector<RooAbsBinning *>::iterator itr = _referenceGrid._grid.begin(); itr != _referenceGrid._grid.end();
201 ++itr) {
202 _referenceGrid._nnuis.push_back((*itr)->numBins() + 1);
203 }
204
205 int nPar = _parList.size();
206 int nDim = _referenceGrid._grid.size();
208 int nRef = _referenceGrid._nref.size();
209 int depth = std::pow(2, nPar);
210
211 if (nPar != nDim) {
212 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPar != nDim"
213 << ": " << nPar << " !=" << nDim << std::endl;
214 assert(0);
215 }
216
217 if (nPdf != nRef) {
218 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPdf != nRef"
219 << ": " << nPdf << " !=" << nRef << std::endl;
220 assert(0);
221 }
222
223 // Transformation matrix for NonLinear settings
224 _M = std::make_unique<TMatrixD>(nPdf, nPdf);
225 _MSqr = std::make_unique<TMatrixD>(depth, depth);
227 TMatrixD M(nPdf, nPdf);
228
230 for (int k = 0; k < nPdf; ++k) {
232 for (int idim = 0; idim < nPar; idim++) {
233 double delta = _referenceGrid._nref[k][idim] - _referenceGrid._nref[0][idim];
234 dm2.push_back(delta);
235 }
236 dm[k] = dm2;
237 }
238
240 for (int idim = 0; idim < nPar; idim++) {
242 xtmp.reserve(_referenceGrid._nnuis[idim]);
243 for (int ix = 0; ix < _referenceGrid._nnuis[idim]; ix++) {
244 xtmp.push_back(ix);
245 }
246 powers.push_back(xtmp);
247 }
248
251 int nCombs = output.size();
252
253 for (int k = 0; k < nPdf; ++k) {
254 int nperm = 0;
255 for (int i = 0; i < nCombs; i++) {
256 double tmpDm = 1.0;
257 for (int ix = 0; ix < nPar; ix++) {
258 double delta = dm[k][ix];
259 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
260 }
261 M(k, nperm) = tmpDm;
262 nperm++;
263 }
264 }
265
266 // M.Print();
267 (*_M) = M.Invert();
268 }
269
270 // Resize transformation vectors
271 _squareVec.resize(std::pow(2, nPar));
272 _squareIdx.resize(std::pow(2, nPar));
273}
274
275//_____________________________________________________________________________
277 : _pdfList(other._pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
278{
279 for (unsigned int i = 0; i < other._grid.size(); i++) {
280 _grid.push_back(other._grid[i]->clone());
281 }
282}
283
284//_____________________________________________________________________________
286{
287 for (RooAbsBinning *binning : _grid) {
288 delete binning;
289 }
290}
291
292//_____________________________________________________________________________
294{
297 thisBoundaries.push_back(bin_x);
298 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
299 _pdfList.add(pdf);
300 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
301 _nref.push_back(thisBoundaryCoordinates);
302}
303
304//_____________________________________________________________________________
306{
309 thisBoundaries.push_back(bin_x);
310 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
311 thisBoundaries.push_back(bin_y);
312 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
313 _pdfList.add(pdf);
314 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
315 _nref.push_back(thisBoundaryCoordinates);
316}
317
318//_____________________________________________________________________________
320{
323 thisBoundaries.push_back(bin_x);
324 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
325 thisBoundaries.push_back(bin_y);
326 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
327 thisBoundaries.push_back(bin_z);
328 thisBoundaryCoordinates.push_back(_grid[2]->array()[bin_z]);
329 _pdfList.add(pdf);
330 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
331 _nref.push_back(thisBoundaryCoordinates);
332}
333
334//_____________________________________________________________________________
336{
338 int nBins = bins.size();
339 thisBoundaryCoordinates.reserve(nBins);
340 for (int i = 0; i < nBins; i++) {
341 thisBoundaryCoordinates.push_back(_grid[i]->array()[bins[i]]);
342 }
343 _pdfList.add(pdf);
344 _pdfMap[bins] = _pdfList.size() - 1;
345 _nref.push_back(thisBoundaryCoordinates);
346}
347
348//_____________________________________________________________________________
350{
351 auto cache = static_cast<CacheElem *>(_cacheMgr.getObj(nullptr, static_cast<RooArgSet const*>(nullptr)));
352 if (cache) {
353 return cache;
354 }
355
356 int nObs = _obsList.size();
358
359 RooAbsReal *null = nullptr;
360 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
361 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
362 vector<RooAbsReal *> myrms(nObs, null);
363 vector<RooAbsReal *> mypos(nObs, null);
364 vector<RooAbsReal *> slope(nPdf * nObs, null);
365 vector<RooAbsReal *> offsets(nPdf * nObs, null);
366 vector<RooAbsReal *> transVar(nPdf * nObs, null);
368
371
372 // fraction parameters
373 RooArgList coefList("coefList"); // fractions multiplied with input pdfs
374 RooArgList coefList2("coefList2"); // fractions multiplied with mean position of observable contribution
375 RooArgList coefList3("coefList3"); // fractions multiplied with rms position of observable contribution
376
377 for (int i = 0; i < 3 * nPdf; ++i) {
378 string fracName = Form("frac_%d", i);
379 double initval = _isPdfMode ? 1.0 : 0.0;
380 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), initval); // to be set later
381
382 fracl.add(*frac);
383 if (i < nPdf) {
384 coefList.add(*static_cast<RooRealVar *>(fracl.at(i)));
385 } else if (i < 2 * nPdf) {
386 coefList2.add(*static_cast<RooRealVar *>(fracl.at(i)));
387 } else {
388 coefList3.add(*static_cast<RooRealVar *>(fracl.at(i)));
389 }
390 ownedComps.add(*static_cast<RooRealVar *>(fracl.at(i)));
391 }
392
393 std::unique_ptr<RooAbsReal> theSum;
394 string sumName = Form("%s_sum", GetName());
395
397 if (_useHorizMorph) {
398 // mean and sigma
400 for (int i = 0; i < nPdf; ++i) {
401 for (int j = 0; j < nObs; ++j) {
402 RooAbsMoment *mom = nObs == 1 ? (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)))
403 : (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)), obsList);
404
405 mom->setLocalNoDirtyInhibit(true);
406 mom->mean()->setLocalNoDirtyInhibit(true);
407
408 sigmarv[sij(i, j)] = mom;
409 meanrv[sij(i, j)] = mom->mean();
410
411 ownedComps.add(*sigmarv[sij(i, j)]);
412 }
413 }
414
415 // slope and offset (to be set later, depend on nuisance parameters)
416 for (int j = 0; j < nObs; ++j) {
417 RooArgList meanList("meanList");
418 RooArgList rmsList("rmsList");
419 for (int i = 0; i < nPdf; ++i) {
420 meanList.add(*meanrv[sij(i, j)]);
421 rmsList.add(*sigmarv[sij(i, j)]);
422 }
423 string myrmsName = Form("%s_rms_%d", GetName(), j);
424 string myposName = Form("%s_pos_%d", GetName(), j);
425 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
426 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
427 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
428 }
429
430 // construction of unit pdfs
431
432 Int_t i = 0;
433 for (auto const *pdf : static_range_cast<Base_t *>(_pdfList)) {
434
435 string pdfName = Form("pdf_%d", i);
436 RooCustomizer cust(*pdf, pdfName.c_str());
437
438 Int_t j = 0;
439 for (auto *var : static_range_cast<RooRealVar *>(obsList)) {
440 // slope and offset formulas
441 string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
442 string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
443
444 slope[sij(i, j)] =
445 new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[sij(i, j)], *myrms[j]));
446 offsets[sij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
447 RooArgList(*meanrv[sij(i, j)], *mypos[j], *slope[sij(i, j)]));
448 ownedComps.add(RooArgSet(*slope[sij(i, j)], *offsets[sij(i, j)]));
449
450 // linear transformations, so pdf can be renormalized easily
451 string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
452 transVar[sij(i, j)] = new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[sij(i, j)],
453 *offsets[sij(i, j)]);
454
455 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
456 // This will prevent the likelihood optimizers from erroneously declaring terms constant
457 transVar[sij(i, j)]->addServerList((RooAbsCollection &)_parList);
458
459 ownedComps.add(*transVar[sij(i, j)]);
460 cust.replaceArg(*var, *transVar[sij(i, j)]);
461 ++j;
462 }
463 transPdf[i] = static_cast<Base_t *>(cust.build());
464 transPdfList.add(*transPdf[i]);
465 ownedComps.add(*transPdf[i]);
466 ++i;
467 }
468 }
469
470 // sum pdf
471 RooArgList const &pdfList = _useHorizMorph ? transPdfList : static_cast<RooArgList const &>(_pdfList);
472 if (_isPdfMode) {
473 theSum = std::make_unique<RooAddPdf>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
474 } else {
475 theSum = std::make_unique<RooRealSumFunc>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
476 }
477
478 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
479 // This will prevent the likelihood optimizers from erroneously declaring terms constant
480 theSum->addServerList((RooAbsCollection &)_parList);
481 theSum->addOwnedComponents(ownedComps);
482
483 // change tracker for fraction parameters
484 std::string trackerName = std::string(GetName()) + "_frac_tracker";
485
486 // Store it in the cache
487 cache = new CacheElem(std::move(theSum),
488 std::make_unique<RooChangeTracker>(trackerName.c_str(), trackerName.c_str(), _parList, true),
489 fracl);
490 _cacheMgr.setObj(nullptr, nullptr, cache, nullptr);
491
492 return cache;
493}
494
496 std::unique_ptr<RooChangeTracker> &&tracker, const RooArgList &flist)
497 : _sum(std::move(sumFunc)), _tracker(std::move(tracker))
498{
499 _frac.add(flist);
500}
501
502//_____________________________________________________________________________
507
508//_____________________________________________________________________________
510
511//_____________________________________________________________________________
513{
514 // Special version of getValV() overrides Base_t::getValV() to save value of current normalization set
515 _curNormSet = set ? const_cast<RooArgSet *>(set) : const_cast<RooArgSet *>(static_cast<RooArgSet const*>(&_obsList));
516 return Base_t::getValV(set);
517}
518
519//_____________________________________________________________________________
521{
522 CacheElem *cache = getCache(nset ? nset : _curNormSet);
523
524 if (cache->_tracker->hasChanged(true)) {
525 cache->calculateFractions(*this, false); // verbose turned off
526 }
527 return cache->_sum.get();
528}
529
530//_____________________________________________________________________________
532{
534
535 if (cache->_tracker->hasChanged(true)) {
536 cache->calculateFractions(*this, false); // verbose turned off
537 }
538
539 double ret = cache->_sum->getVal(_obsList.nset());
540
541 return ret;
542}
543
544//_____________________________________________________________________________
546{
547 return static_cast<RooRealVar *>(_frac.at(i));
548}
549
550//_____________________________________________________________________________
552{
553 return static_cast<RooRealVar *>(_frac.at(i));
554}
555
556//_____________________________________________________________________________
558{
559 int nPdf = self._pdfList.size();
560 int nPar = self._parList.size();
561
562 double fracLinear(1.);
563 double fracNonLinear(1.);
564
565 if (self._setting == NonLinear || self._setting == NonLinearLinFractions || self._setting == NonLinearPosFractions) {
566 // Calculate the delta vector
568 for (int idim = 0; idim < nPar; idim++) {
569 double delta = (static_cast<RooRealVar *>(self._parList.at(idim)))->getVal() - self._referenceGrid._nref[0][idim];
570 dm2.push_back(delta);
571 }
572
574 for (int idim = 0; idim < nPar; idim++) {
576 xtmp.reserve(self._referenceGrid._nnuis[idim]);
577 for (int ix = 0; ix < self._referenceGrid._nnuis[idim]; ix++) {
578 xtmp.push_back(ix);
579 }
580 powers.push_back(xtmp);
581 }
582
585 int nCombs = output.size();
586
588
589 int nperm = 0;
590 for (int i = 0; i < nCombs; i++) {
591 double tmpDm = 1.0;
592 for (int ix = 0; ix < nPar; ix++) {
593 double delta = dm2[ix];
594 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
595 }
597 nperm++;
598 }
599
600 double sumposfrac = 0.0;
601 for (int i = 0; i < nPdf; ++i) {
602 double ffrac = 0.0;
603
604 for (int j = 0; j < nPdf; ++j) {
605 ffrac += (*self._M)(j, i) * deltavec[j] * fracNonLinear;
606 }
607
608 if (ffrac >= 0) {
609 sumposfrac += ffrac;
610 }
611
612 // fractions for pdf
613 if (self._setting != NonLinearLinFractions) {
614 const_cast<RooRealVar *>(frac(i))->setVal(ffrac);
615 }
616
617 // fractions for rms and mean
618 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(ffrac); // need to add up
619 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(ffrac); // need to add up
620
621 if (verbose) {
622 std::cout << "NonLinear fraction " << ffrac << std::endl;
623 frac(i)->Print();
624 frac(nPdf + i)->Print();
625 frac(2 * nPdf + i)->Print();
626 }
627 }
628
629 if (self._setting == NonLinearPosFractions) {
630 for (int i = 0; i < nPdf; ++i) {
631 if (frac(i)->getVal() < 0)
632 const_cast<RooRealVar *>(frac(i))->setVal(0.);
633 const_cast<RooRealVar *>(frac(i))->setVal(frac(i)->getVal() / sumposfrac);
634 }
635 }
636 }
637
638 if (self._setting == Linear || self._setting == NonLinearLinFractions) {
639 // zero all fractions
640 // for (int i = 0; i < 3*nPdf; ++i) {
641 for (int i = 0; i < nPdf; ++i) {
642 double initval = 0;
643 const_cast<RooRealVar *>(frac(i))->setVal(initval);
644 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(initval);
645 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(initval);
646 }
647
648 std::vector<double> mtmp;
649
650 // loop over parList
651 for (auto *m : static_range_cast<RooRealVar *>(self._parList)) {
652 mtmp.push_back(m->getVal());
653 }
654
655 self.findShape(mtmp); // this sets _squareVec and _squareIdx quantities
656
657 int depth = std::pow(2, nPar);
659
660 int nperm = 0;
661
663 xtmp.reserve(nPar);
664 for (int ix = 0; ix < nPar; ix++) {
665 xtmp.push_back(ix);
666 }
667
668 for (int iperm = 1; iperm <= nPar; ++iperm) {
669 do {
670 double dtmp = mtmp[xtmp[0]] - self._squareVec[0][xtmp[0]];
671 for (int itmp = 1; itmp < iperm; ++itmp) {
672 dtmp *= mtmp[xtmp[itmp]] - self._squareVec[0][xtmp[itmp]];
673 }
674 deltavec[nperm + 1] = dtmp;
675 nperm++;
677 }
678
679 double origFrac1(0.);
680 double origFrac2(0.);
681 for (int i = 0; i < depth; ++i) {
682 double ffrac = 0.;
683 for (int j = 0; j < depth; ++j) {
684 ffrac += (*self._MSqr)(j, i) * deltavec[j] * fracLinear;
685 }
686
687 // set fractions for pdf
688 origFrac1 = frac(self._squareIdx[i])->getVal(); // already set in case of smoothlinear
689 const_cast<RooRealVar *>(frac(self._squareIdx[i]))->setVal(origFrac1 + ffrac); // need to add up
690
691 // set fractions for rms and mean
692 if (self._setting != NonLinearLinFractions) {
693 origFrac2 =
694 frac(nPdf + self._squareIdx[i])->getVal(); // already set in case of smoothlinear
695 const_cast<RooRealVar *>(frac(nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
696 const_cast<RooRealVar *>(frac(2 * nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
697 }
698
699 if (verbose) {
700 std::cout << "Linear fraction " << ffrac << std::endl;
701 frac(self._squareIdx[i])->Print();
702 frac(nPdf + self._squareIdx[i])->Print();
703 frac(2 * nPdf + self._squareIdx[i])->Print();
704 }
705 }
706 }
707}
708
709//_____________________________________________________________________________
711{
712 int nPar = _parList.size();
713 int nRef = _referenceGrid._nref.size();
714
715 // Find hypercube enclosing the location to morph to
716 // bool isEnclosed = true;
717 // for (int i = 0; i < nPar; i++) {
718 // if (x[i] < _referenceGrid._grid[i]->lowBound())
719 // isEnclosed = false;
720 // if (x[i] > _referenceGrid._grid[i]->highBound())
721 // isEnclosed = false;
722 // }
723
724 // std::cout << "isEnclosed = " << isEnclosed << std::endl;
725
726 int depth = std::pow(2, nPar);
727
728 vector<vector<double>> boundaries(nPar);
729 for (int idim = 0; idim < nPar; idim++) {
730 int bin = _referenceGrid._grid[idim]->binNumber(x[idim]);
731 double lo = _referenceGrid._grid[idim]->binLow(bin);
732 double hi = _referenceGrid._grid[idim]->binHigh(bin);
733 boundaries[idim].push_back(lo);
734 boundaries[idim].push_back(hi);
735 }
736
740
741 for (int isq = 0; isq < depth; isq++) {
742 for (int iref = 0; iref < nRef; iref++) {
745 break;
746 }
747 }
748 }
749
750 // std::cout << std::endl;
751
752 // for (int isq = 0; isq < _squareVec.size(); isq++) {
753 // std::cout << _squareIdx[isq];
754 // std::cout << " (";
755 // for (int isqq = 0; isqq < _squareVec[isq].size(); isqq++) {
756 // std::cout << _squareVec[isq][isqq] << ((isqq<_squareVec[isq].size()-1)?",":"");
757 // }
758 // std::cout << ") ";
759 // }
760
761 // construct transformation matrix for linear extrapolation
762 TMatrixD M(depth, depth);
763
765 xtmp.reserve(nPar);
766 for (int ix = 0; ix < nPar; ix++) {
767 xtmp.push_back(ix);
768 }
769
770 for (int k = 0; k < depth; ++k) {
771 M(k, 0) = 1.0;
772
773 int nperm = 0;
775
776 for (int iperm = 1; iperm <= nPar; ++iperm) {
777 do {
778 double dtmp = _squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
779 for (int itmp = 1; itmp < iperm; ++itmp) {
781 }
782 M(k, nperm + 1) = dtmp;
783 nperm++;
785 }
786 }
787
788 // M.Print();
789 (*_MSqr) = M.Invert();
790}
791
792//_____________________________________________________________________________
794{
795 if (allVars.size() == 1) {
796 RooAbsReal *temp = const_cast<RooMomentMorphFuncND *>(this);
797 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator");
798 int nbins = (static_cast<RooRealVar *>(allVars.first()))->numBins();
799 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins", nbins);
800 return true;
801 } else {
802 std::cout << "Currently BinIntegrator only knows how to deal with 1-d " << std::endl;
803 return false;
804 }
805 return false;
806}
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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:2489
const_iterator begin() const
const_iterator end() const
Abstract base class for RooRealVar binning definitions.
Abstract container object that can hold multiple RooAbsArg objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * first() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h: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.
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
Calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition RooAddition.h:27
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Implements a RooAbsBinning in terms of an array of boundary values, posing no constraints on the choi...
Definition RooBinning.h:27
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
bool setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Represents a constant real-valued object.
Definition RooConstVar.h:23
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooLinearVar is the most general form of a derived real-valued object that can be used by RooRealInte...
void calculateFractions(const RooMomentMorphFuncND &self, bool verbose=true) const
std::unique_ptr< RooChangeTracker > _tracker
std::unique_ptr< RooAbsReal > _sum
RooArgList containedArgs(Action) override
CacheElem(std::unique_ptr< RooAbsReal > &&sumFunc, std::unique_ptr< RooChangeTracker > &&tracker, const RooArgList &flist)
void addBinning(const RooAbsBinning &binning)
std::vector< RooAbsBinning * > _grid
std::vector< std::vector< double > > _nref
void addPdf(const RooAbsReal &func, int bin_x)
RooObjCacheManager _cacheMgr
! Transient cache manager
std::unique_ptr< TMatrixD > _MSqr
void findShape(const std::vector< double > &x) const
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsReal * sumFunc(const RooArgSet *nset)
CacheElem * getCache(const RooArgSet *nset) const
std::unique_ptr< TMatrixD > _M
RooArgSet * _curNormSet
! Transient cache manager
bool setBinIntegrator(RooArgSet &allVars)
double getValV(const RooArgSet *set=nullptr) const override
Return value of object.
std::vector< std::vector< double > > _squareVec
int sij(const int &i, const int &j) const
std::vector< int > _squareIdx
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooCategory & method1D()
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Double_t x[n]
Definition legend1.C:17
void cartesianProduct(std::vector< std::vector< T > > &out, std::vector< std::vector< T > > &in)
Definition Algorithms.h:22
bool nextCombination(const Iterator first, Iterator k, const Iterator last)
Definition Algorithms.h:64
TMarker m
Definition textangle.C:8
static void output()