Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorphFuncND.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2023, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11/** \class RooMomentMorphFuncND
12 \ingroup Roofit
13
14**/
15
17
18#include <RooAbsCategory.h>
19#include <RooAbsMoment.h>
20#include <RooAddPdf.h>
21#include <RooAddition.h>
22#include <RooChangeTracker.h>
23#include <RooConstVar.h>
24#include <RooCustomizer.h>
25#include <RooFormulaVar.h>
26#include <RooLinearVar.h>
27#include <RooMoment.h>
28#include <RooNumIntConfig.h>
29#include <RooRealSumFunc.h>
30#include <RooRealVar.h>
31
33
34#include <Riostream.h>
35
36#include <TMap.h>
37#include <TMath.h>
38#include <TVector.h>
39
40#include <algorithm>
41#include <map>
42
43using std::cout, std::endl, std::string, std::vector;
44
47
48//_____________________________________________________________________________
49RooMomentMorphFuncND::RooMomentMorphFuncND() : _cacheMgr(this, 10, true, true), _setting(RooMomentMorphFuncND::Linear), _useHorizMorph(true)
50{
52}
53
54//_____________________________________________________________________________
55RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, const RooArgList &parList, const RooArgList &obsList,
56 const Grid2 &referenceGrid, const Setting &setting)
58 _cacheMgr(this, 10, true, true),
59 _parList("parList", "List of morph parameters", this),
60 _obsList("obsList", "List of observables", this),
61 _referenceGrid(referenceGrid),
62 _pdfList("pdfList", "List of pdfs", this),
63 _setting(setting),
64 _useHorizMorph(true)
65{
66 // morph parameters
68
69 // observables
71
73
74 // general initialization
75 initialize();
76
78}
79
80//_____________________________________________________________________________
81RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
82 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
84 _cacheMgr(this, 10, true, true),
85 _parList("parList", "List of morph parameters", this),
86 _obsList("obsList", "List of observables", this),
87 _pdfList("pdfList", "List of pdfs", this),
88 _setting(setting),
89 _useHorizMorph(true)
90{
91 // make reference grid
92 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
94
95 for (int i = 0; i < mrefpoints.GetNrows(); ++i) {
96 for (int j = 0; j < grid.numBoundaries(); ++j) {
97 if (mrefpoints[i] == grid.array()[j]) {
98 _referenceGrid.addPdf(*static_cast<Base_t *>(pdfList.at(i)), j);
99 break;
100 }
101 }
102 }
103
105
106 // morph parameters
107 RooArgList parList;
108 parList.add(_m);
109 _parList.addTyped<RooAbsReal>(parList);
110
111 // observables
112 _obsList.addTyped<RooAbsReal>(varList);
113
114 // general initialization
115 initialize();
116
118}
119
120//_____________________________________________________________________________
121RooMomentMorphFuncND::RooMomentMorphFuncND(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
122 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
124 _cacheMgr(this, 10, true, true),
125 _parList("parList", "List of morph parameters", this),
126 _obsList("obsList", "List of observables", this),
127 _pdfList("pdfList", "List of pdfs", this),
128 _setting(setting),
129 _useHorizMorph(true)
130{
131 // make reference grid
132 TVectorD mrefpoints(mrefList.size());
133 Int_t i = 0;
134 for (auto *mref : mrefList) {
135 if (!dynamic_cast<RooAbsReal *>(mref)) {
136 coutE(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
137 << " is not of type RooAbsReal" << endl;
138 throw string("RooMomentMorphFuncND::ctor() ERROR mref is not of type RooAbsReal");
139 }
140 if (!dynamic_cast<RooConstVar *>(mref)) {
141 coutW(InputArguments) << "RooMomentMorphFuncND::ctor(" << GetName() << ") WARNING mref point " << i
142 << " is not a constant, taking a snapshot of its value" << endl;
143 }
144 mrefpoints[i] = static_cast<RooAbsReal *>(mref)->getVal();
145 i++;
146 }
147
148 RooBinning grid(mrefpoints.GetNrows() - 1, mrefpoints.GetMatrixArray());
150 for (i = 0; i < mrefpoints.GetNrows(); ++i) {
151 for (int j = 0; j < grid.numBoundaries(); ++j) {
152 if (mrefpoints[i] == grid.array()[j]) {
153 _referenceGrid.addPdf(static_cast<Base_t &>(pdfList[i]), j);
154 break;
155 }
156 }
157 }
158
160
161 // morph parameters
162 RooArgList parList;
163 parList.add(_m);
164 _parList.addTyped<RooAbsReal>(parList);
165
166 // observables
167 _obsList.addTyped<RooAbsReal>(varList);
168
169 // general initialization
170 initialize();
171
173}
174
175//_____________________________________________________________________________
178 _cacheMgr(other._cacheMgr, this),
179 _parList("parList", this, other._parList),
180 _obsList("obsList", this, other._obsList),
181 _referenceGrid(other._referenceGrid),
182 _pdfList("pdfList", this, other._pdfList),
183 _setting(other._setting),
184 _useHorizMorph(other._useHorizMorph),
185 _isPdfMode{other._isPdfMode}
186{
187 // general initialization
188 initialize();
189
191}
192
193//_____________________________________________________________________________
195{
197}
198
199//_____________________________________________________________________________
201{
202 for (vector<RooAbsBinning *>::iterator itr = _referenceGrid._grid.begin(); itr != _referenceGrid._grid.end();
203 ++itr) {
204 _referenceGrid._nnuis.push_back((*itr)->numBins() + 1);
205 }
206
207 int nPar = _parList.size();
208 int nDim = _referenceGrid._grid.size();
209 int nPdf = _referenceGrid._pdfList.size();
210 int nRef = _referenceGrid._nref.size();
211 int depth = std::pow(2, nPar);
212
213 if (nPar != nDim) {
214 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPar != nDim"
215 << ": " << nPar << " !=" << nDim << endl;
216 assert(0);
217 }
218
219 if (nPdf != nRef) {
220 coutE(InputArguments) << "RooMomentMorphFuncND::initialize(" << GetName() << ") ERROR: nPdf != nRef"
221 << ": " << nPdf << " !=" << nRef << endl;
222 assert(0);
223 }
224
225 // Transformation matrix for NonLinear settings
226 _M = std::make_unique<TMatrixD>(nPdf, nPdf);
227 _MSqr = std::make_unique<TMatrixD>(depth, depth);
229 TMatrixD M(nPdf, nPdf);
230
231 vector<vector<double>> dm(nPdf);
232 for (int k = 0; k < nPdf; ++k) {
233 vector<double> dm2;
234 for (int idim = 0; idim < nPar; idim++) {
235 double delta = _referenceGrid._nref[k][idim] - _referenceGrid._nref[0][idim];
236 dm2.push_back(delta);
237 }
238 dm[k] = dm2;
239 }
240
241 vector<vector<int>> powers;
242 for (int idim = 0; idim < nPar; idim++) {
243 vector<int> xtmp;
244 xtmp.reserve(_referenceGrid._nnuis[idim]);
245 for (int ix = 0; ix < _referenceGrid._nnuis[idim]; ix++) {
246 xtmp.push_back(ix);
247 }
248 powers.push_back(xtmp);
249 }
250
251 vector<vector<int>> output;
253 int nCombs = output.size();
254
255 for (int k = 0; k < nPdf; ++k) {
256 int nperm = 0;
257 for (int i = 0; i < nCombs; i++) {
258 double tmpDm = 1.0;
259 for (int ix = 0; ix < nPar; ix++) {
260 double delta = dm[k][ix];
261 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
262 }
263 M(k, nperm) = tmpDm;
264 nperm++;
265 }
266 }
267
268 // M.Print();
269 (*_M) = M.Invert();
270 }
271
272 // Resize transformation vectors
273 _squareVec.resize(std::pow(2, nPar));
274 _squareIdx.resize(std::pow(2, nPar));
275}
276
277//_____________________________________________________________________________
279 : _pdfList(other._pdfList), _pdfMap(other._pdfMap), _nref(other._nref)
280{
281 for (unsigned int i = 0; i < other._grid.size(); i++) {
282 _grid.push_back(other._grid[i]->clone());
283 }
284}
285
286//_____________________________________________________________________________
288{
289 for (RooAbsBinning *binning : _grid) {
290 delete binning;
291 }
292}
293
294//_____________________________________________________________________________
296{
297 vector<int> thisBoundaries;
298 vector<double> thisBoundaryCoordinates;
299 thisBoundaries.push_back(bin_x);
300 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
301 _pdfList.add(pdf);
302 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
303 _nref.push_back(thisBoundaryCoordinates);
304}
305
306//_____________________________________________________________________________
308{
309 vector<int> thisBoundaries;
310 vector<double> thisBoundaryCoordinates;
311 thisBoundaries.push_back(bin_x);
312 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
313 thisBoundaries.push_back(bin_y);
314 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
315 _pdfList.add(pdf);
316 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
317 _nref.push_back(thisBoundaryCoordinates);
318}
319
320//_____________________________________________________________________________
321void RooMomentMorphFuncND::Grid2::addPdf(const RooMomentMorphFuncND::Base_t &pdf, int bin_x, int bin_y, int bin_z)
322{
323 vector<int> thisBoundaries;
324 vector<double> thisBoundaryCoordinates;
325 thisBoundaries.push_back(bin_x);
326 thisBoundaryCoordinates.push_back(_grid[0]->array()[bin_x]);
327 thisBoundaries.push_back(bin_y);
328 thisBoundaryCoordinates.push_back(_grid[1]->array()[bin_y]);
329 thisBoundaries.push_back(bin_z);
330 thisBoundaryCoordinates.push_back(_grid[2]->array()[bin_z]);
331 _pdfList.add(pdf);
332 _pdfMap[thisBoundaries] = _pdfList.size() - 1;
333 _nref.push_back(thisBoundaryCoordinates);
334}
335
336//_____________________________________________________________________________
338{
339 vector<double> thisBoundaryCoordinates;
340 int nBins = bins.size();
341 thisBoundaryCoordinates.reserve(nBins);
342 for (int i = 0; i < nBins; i++) {
343 thisBoundaryCoordinates.push_back(_grid[i]->array()[bins[i]]);
344 }
345 _pdfList.add(pdf);
346 _pdfMap[bins] = _pdfList.size() - 1;
347 _nref.push_back(thisBoundaryCoordinates);
348}
349
350//_____________________________________________________________________________
352{
353 auto cache = static_cast<CacheElem *>(_cacheMgr.getObj(nullptr, static_cast<RooArgSet const*>(nullptr)));
354 if (cache) {
355 return cache;
356 }
357
358 int nObs = _obsList.size();
359 int nPdf = _referenceGrid._pdfList.size();
360
361 RooAbsReal *null = nullptr;
362 vector<RooAbsReal *> meanrv(nPdf * nObs, null);
363 vector<RooAbsReal *> sigmarv(nPdf * nObs, null);
364 vector<RooAbsReal *> myrms(nObs, null);
365 vector<RooAbsReal *> mypos(nObs, null);
366 vector<RooAbsReal *> slope(nPdf * nObs, null);
367 vector<RooAbsReal *> offsets(nPdf * nObs, null);
368 vector<RooAbsReal *> transVar(nPdf * nObs, null);
369 vector<RooAbsReal *> transPdf(nPdf, null);
370
371 RooArgSet ownedComps;
372 RooArgList fracl;
373
374 // fraction parameters
375 RooArgList coefList("coefList"); // fractions multiplied with input pdfs
376 RooArgList coefList2("coefList2"); // fractions multiplied with mean position of observable contribution
377 RooArgList coefList3("coefList3"); // fractions multiplied with rms position of observable contribution
378
379 for (int i = 0; i < 3 * nPdf; ++i) {
380 string fracName = Form("frac_%d", i);
381 double initval = _isPdfMode ? 1.0 : 0.0;
382 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), initval); // to be set later
383
384 fracl.add(*frac);
385 if (i < nPdf) {
386 coefList.add(*static_cast<RooRealVar *>(fracl.at(i)));
387 } else if (i < 2 * nPdf) {
388 coefList2.add(*static_cast<RooRealVar *>(fracl.at(i)));
389 } else {
390 coefList3.add(*static_cast<RooRealVar *>(fracl.at(i)));
391 }
392 ownedComps.add(*static_cast<RooRealVar *>(fracl.at(i)));
393 }
394
395 std::unique_ptr<RooAbsReal> theSum;
396 string sumName = Form("%s_sum", GetName());
397
398 RooArgList transPdfList;
399 if (_useHorizMorph) {
400 // mean and sigma
401 RooArgList obsList(_obsList);
402 for (int i = 0; i < nPdf; ++i) {
403 for (int j = 0; j < nObs; ++j) {
404 RooAbsMoment *mom = nObs == 1 ? (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)))
405 : (static_cast<Base_t *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*obsList.at(j)), obsList);
406
407 mom->setLocalNoDirtyInhibit(true);
408 mom->mean()->setLocalNoDirtyInhibit(true);
409
410 sigmarv[sij(i, j)] = mom;
411 meanrv[sij(i, j)] = mom->mean();
412
413 ownedComps.add(*sigmarv[sij(i, j)]);
414 }
415 }
416
417 // slope and offset (to be set later, depend on nuisance parameters)
418 for (int j = 0; j < nObs; ++j) {
419 RooArgList meanList("meanList");
420 RooArgList rmsList("rmsList");
421 for (int i = 0; i < nPdf; ++i) {
422 meanList.add(*meanrv[sij(i, j)]);
423 rmsList.add(*sigmarv[sij(i, j)]);
424 }
425 string myrmsName = Form("%s_rms_%d", GetName(), j);
426 string myposName = Form("%s_pos_%d", GetName(), j);
427 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
428 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList3);
429 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
430 }
431
432 // construction of unit pdfs
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
472 // sum pdf
473 RooArgList const &pdfList = _useHorizMorph ? transPdfList : static_cast<RooArgList const &>(_pdfList);
474 if (_isPdfMode) {
475 theSum = std::make_unique<RooAddPdf>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
476 } else {
477 theSum = std::make_unique<RooRealSumFunc>(sumName.c_str(), sumName.c_str(), pdfList, coefList);
478 }
479
480 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
481 // This will prevent the likelihood optimizers from erroneously declaring terms constant
482 theSum->addServerList((RooAbsCollection &)_parList);
483 theSum->addOwnedComponents(ownedComps);
484
485 // change tracker for fraction parameters
486 std::string trackerName = std::string(GetName()) + "_frac_tracker";
487
488 // Store it in the cache
489 cache = new CacheElem(std::move(theSum),
490 std::make_unique<RooChangeTracker>(trackerName.c_str(), trackerName.c_str(), _parList, true),
491 fracl);
492 _cacheMgr.setObj(nullptr, nullptr, cache, nullptr);
493
494 return cache;
495}
496
498 std::unique_ptr<RooChangeTracker> &&tracker, const RooArgList &flist)
499 : _sum(std::move(sumFunc)), _tracker(std::move(tracker))
500{
501 _frac.add(flist);
502}
503
504//_____________________________________________________________________________
506{
507 return RooArgList(*_sum, *_tracker);
508}
509
510//_____________________________________________________________________________
512
513//_____________________________________________________________________________
515{
516 // Special version of getValV() overrides Base_t::getValV() to save value of current normalization set
517 _curNormSet = set ? const_cast<RooArgSet *>(set) : const_cast<RooArgSet *>(static_cast<RooArgSet const*>(&_obsList));
518 return Base_t::getValV(set);
519}
520
521//_____________________________________________________________________________
523{
524 CacheElem *cache = getCache(nset ? nset : _curNormSet);
525
526 if (cache->_tracker->hasChanged(true)) {
527 cache->calculateFractions(*this, false); // verbose turned off
528 }
529 return cache->_sum.get();
530}
531
532//_____________________________________________________________________________
534{
536
537 if (cache->_tracker->hasChanged(true)) {
538 cache->calculateFractions(*this, false); // verbose turned off
539 }
540
541 double ret = cache->_sum->getVal(_obsList.nset());
542
543 return ret;
544}
545
546//_____________________________________________________________________________
548{
549 return static_cast<RooRealVar *>(_frac.at(i));
550}
551
552//_____________________________________________________________________________
554{
555 return static_cast<RooRealVar *>(_frac.at(i));
556}
557
558//_____________________________________________________________________________
560{
561 int nPdf = self._pdfList.size();
562 int nPar = self._parList.size();
563
564 double fracLinear(1.);
565 double fracNonLinear(1.);
566
568 // Calculate the delta vector
569 vector<double> dm2;
570 for (int idim = 0; idim < nPar; idim++) {
571 double delta = (static_cast<RooRealVar *>(self._parList.at(idim)))->getVal() - self._referenceGrid._nref[0][idim];
572 dm2.push_back(delta);
573 }
574
575 vector<vector<int>> powers;
576 for (int idim = 0; idim < nPar; idim++) {
577 vector<int> xtmp;
578 xtmp.reserve(self._referenceGrid._nnuis[idim]);
579 for (int ix = 0; ix < self._referenceGrid._nnuis[idim]; ix++) {
580 xtmp.push_back(ix);
581 }
582 powers.push_back(xtmp);
583 }
584
585 vector<vector<int>> output;
587 int nCombs = output.size();
588
589 vector<double> deltavec(nPdf, 1.0);
590
591 int nperm = 0;
592 for (int i = 0; i < nCombs; i++) {
593 double tmpDm = 1.0;
594 for (int ix = 0; ix < nPar; ix++) {
595 double delta = dm2[ix];
596 tmpDm *= std::pow(delta, static_cast<double>(output[i][ix]));
597 }
598 deltavec[nperm] = tmpDm;
599 nperm++;
600 }
601
602 double sumposfrac = 0.0;
603 for (int i = 0; i < nPdf; ++i) {
604 double ffrac = 0.0;
605
606 for (int j = 0; j < nPdf; ++j) {
607 ffrac += (*self._M)(j, i) * deltavec[j] * fracNonLinear;
608 }
609
610 if (ffrac >= 0) {
611 sumposfrac += ffrac;
612 }
613
614 // fractions for pdf
615 if (self._setting != NonLinearLinFractions) {
616 const_cast<RooRealVar *>(frac(i))->setVal(ffrac);
617 }
618
619 // fractions for rms and mean
620 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(ffrac); // need to add up
621 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(ffrac); // need to add up
622
623 if (verbose) {
624 cout << "NonLinear fraction " << ffrac << endl;
625 frac(i)->Print();
626 frac(nPdf + i)->Print();
627 frac(2 * nPdf + i)->Print();
628 }
629 }
630
631 if (self._setting == NonLinearPosFractions) {
632 for (int i = 0; i < nPdf; ++i) {
633 if (frac(i)->getVal() < 0)
634 const_cast<RooRealVar *>(frac(i))->setVal(0.);
635 const_cast<RooRealVar *>(frac(i))->setVal(frac(i)->getVal() / sumposfrac);
636 }
637 }
638 }
639
640 if (self._setting == Linear || self._setting == NonLinearLinFractions) {
641 // zero all fractions
642 // for (int i = 0; i < 3*nPdf; ++i) {
643 for (int i = 0; i < nPdf; ++i) {
644 double initval = 0;
645 const_cast<RooRealVar *>(frac(i))->setVal(initval);
646 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(initval);
647 const_cast<RooRealVar *>(frac(2 * nPdf + i))->setVal(initval);
648 }
649
650 std::vector<double> mtmp;
651
652 // loop over parList
653 for (auto *m : static_range_cast<RooRealVar *>(self._parList)) {
654 mtmp.push_back(m->getVal());
655 }
656
657 self.findShape(mtmp); // this sets _squareVec and _squareIdx quantities
658
659 int depth = std::pow(2, nPar);
660 vector<double> deltavec(depth, 1.0);
661
662 int nperm = 0;
663
664 vector<int> xtmp;
665 xtmp.reserve(nPar);
666 for (int ix = 0; ix < nPar; ix++) {
667 xtmp.push_back(ix);
668 }
669
670 for (int iperm = 1; iperm <= nPar; ++iperm) {
671 do {
672 double dtmp = mtmp[xtmp[0]] - self._squareVec[0][xtmp[0]];
673 for (int itmp = 1; itmp < iperm; ++itmp) {
674 dtmp *= mtmp[xtmp[itmp]] - self._squareVec[0][xtmp[itmp]];
675 }
676 deltavec[nperm + 1] = dtmp;
677 nperm++;
678 } while (RooFit::Detail::nextCombination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
679 }
680
681 double origFrac1(0.);
682 double origFrac2(0.);
683 for (int i = 0; i < depth; ++i) {
684 double ffrac = 0.;
685 for (int j = 0; j < depth; ++j) {
686 ffrac += (*self._MSqr)(j, i) * deltavec[j] * fracLinear;
687 }
688
689 // set fractions for pdf
690 origFrac1 = frac(self._squareIdx[i])->getVal(); // already set in case of smoothlinear
691 const_cast<RooRealVar *>(frac(self._squareIdx[i]))->setVal(origFrac1 + ffrac); // need to add up
692
693 // set fractions for rms and mean
694 if (self._setting != NonLinearLinFractions) {
695 origFrac2 =
696 frac(nPdf + self._squareIdx[i])->getVal(); // already set in case of smoothlinear
697 const_cast<RooRealVar *>(frac(nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
698 const_cast<RooRealVar *>(frac(2 * nPdf + self._squareIdx[i]))->setVal(origFrac2 + ffrac); // need to add up
699 }
700
701 if (verbose) {
702 cout << "Linear fraction " << ffrac << endl;
703 frac(self._squareIdx[i])->Print();
704 frac(nPdf + self._squareIdx[i])->Print();
705 frac(2 * nPdf + self._squareIdx[i])->Print();
706 }
707 }
708 }
709}
710
711//_____________________________________________________________________________
712void RooMomentMorphFuncND::findShape(const vector<double> &x) const
713{
714 int nPar = _parList.size();
715 int nRef = _referenceGrid._nref.size();
716
717 // Find hypercube enclosing the location to morph to
718 // bool isEnclosed = true;
719 // for (int i = 0; i < nPar; i++) {
720 // if (x[i] < _referenceGrid._grid[i]->lowBound())
721 // isEnclosed = false;
722 // if (x[i] > _referenceGrid._grid[i]->highBound())
723 // isEnclosed = false;
724 // }
725
726 // cout << "isEnclosed = " << isEnclosed << endl;
727
728 int depth = std::pow(2, nPar);
729
730 vector<vector<double>> boundaries(nPar);
731 for (int idim = 0; idim < nPar; idim++) {
732 int bin = _referenceGrid._grid[idim]->binNumber(x[idim]);
733 double lo = _referenceGrid._grid[idim]->binLow(bin);
734 double hi = _referenceGrid._grid[idim]->binHigh(bin);
735 boundaries[idim].push_back(lo);
736 boundaries[idim].push_back(hi);
737 }
738
739 vector<vector<double>> output;
742
743 for (int isq = 0; isq < depth; isq++) {
744 for (int iref = 0; iref < nRef; iref++) {
745 if (_squareVec[isq] == _referenceGrid._nref[iref]) {
746 _squareIdx[isq] = iref;
747 break;
748 }
749 }
750 }
751
752 // cout << endl;
753
754 // for (int isq = 0; isq < _squareVec.size(); isq++) {
755 // cout << _squareIdx[isq];
756 // cout << " (";
757 // for (int isqq = 0; isqq < _squareVec[isq].size(); isqq++) {
758 // cout << _squareVec[isq][isqq] << ((isqq<_squareVec[isq].size()-1)?",":"");
759 // }
760 // cout << ") ";
761 // }
762
763 // construct transformation matrix for linear extrapolation
764 TMatrixD M(depth, depth);
765
766 vector<int> xtmp;
767 xtmp.reserve(nPar);
768 for (int ix = 0; ix < nPar; ix++) {
769 xtmp.push_back(ix);
770 }
771
772 for (int k = 0; k < depth; ++k) {
773 M(k, 0) = 1.0;
774
775 int nperm = 0;
776 vector<double> squareBase = _squareVec[0];
777
778 for (int iperm = 1; iperm <= nPar; ++iperm) {
779 do {
780 double dtmp = _squareVec[k][xtmp[0]] - squareBase[xtmp[0]];
781 for (int itmp = 1; itmp < iperm; ++itmp) {
782 dtmp *= _squareVec[k][xtmp[itmp]] - squareBase[xtmp[itmp]];
783 }
784 M(k, nperm + 1) = dtmp;
785 nperm++;
786 } while (RooFit::Detail::nextCombination(xtmp.begin(), xtmp.begin() + iperm, xtmp.end()));
787 }
788 }
789
790 // M.Print();
791 (*_MSqr) = M.Invert();
792}
793
794//_____________________________________________________________________________
796{
797 if (allVars.size() == 1) {
798 RooAbsReal *temp = const_cast<RooMomentMorphFuncND *>(this);
799 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator");
800 int nbins = (static_cast<RooRealVar *>(allVars.first()))->numBins();
801 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins", nbins);
802 return true;
803 } else {
804 cout << "Currently BinIntegrator only knows how to deal with 1-d " << endl;
805 return false;
806 }
807 return false;
808}
#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:2489
void setLocalNoDirtyInhibit(bool flag) const
Definition RooAbsArg.h:700
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.
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...
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.
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
RooAbsMoment * sigma(RooRealVar &obs)
Definition RooAbsReal.h:368
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
Implements a RooAbsBinning in terms of an array of boundary values, posing no constraints on the choi...
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.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Represents a constant real-valued object.
Definition RooConstVar.h:23
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void 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
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
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:73
Element * GetMatrixArray()
Definition TVectorT.h:76
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()