Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooMomentMorphFunc.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 RooMomentMorphFunc
12 \ingroup Roofit
13
14**/
15
16#include "Riostream.h"
17
18#include "RooMomentMorphFunc.h"
19#include "RooAbsCategory.h"
20#include "RooRealConstant.h"
21#include "RooRealVar.h"
22#include "RooFormulaVar.h"
23#include "RooCustomizer.h"
24#include "RooRealSumFunc.h"
25#include "RooAddition.h"
26#include "RooMoment.h"
27#include "RooLinearVar.h"
28#include "RooChangeTracker.h"
29
30#include "TMath.h"
31#include "TH1.h"
32
33using std::string, std::vector;
34
35
36//_____________________________________________________________________________
38 : _cacheMgr(this, 10, true, true)
39{
40}
41
42//_____________________________________________________________________________
44 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
45 : RooAbsReal(name, title),
46 _cacheMgr(this, 10, true, true),
47 m("m", "m", this, _m),
48 _varList("varList", "List of variables", this),
49 _pdfList("pdfList", "List of pdfs", this),
50 _mref(new TVectorD(mrefpoints)),
51 _setting(setting)
52{
53 // observables
55
56 // reference p.d.f.s
57 _pdfList.addTyped<RooAbsPdf>(pdfList);
58
59 // initialization
60 initialize();
61}
62
63//_____________________________________________________________________________
65 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
66 : RooAbsReal(name, title),
67 _cacheMgr(this, 10, true, true),
68 m("m", "m", this, _m),
69 _varList("varList", "List of variables", this),
70 _pdfList("pdfList", "List of pdfs", this),
71 _mref(new TVectorD(mrefList.size())),
72 _setting(setting)
73{
74 // observables
76
77 // reference p.d.f.s
78 _pdfList.addTyped<RooAbsPdf>(pdfList);
79
80 // reference points in m
81
82 Int_t i = 0;
83 for (auto *mref : mrefList) {
84 if (!dynamic_cast<RooAbsReal *>(mref)) {
85 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
86 << " is not of type RooAbsReal" << std::endl;
87 throw string("RooPolyMorh::ctor() ERROR mref is not of type RooAbsReal");
88 }
89 if (!dynamic_cast<RooConstVar *>(mref)) {
90 coutW(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") WARNING mref point " << i
91 << " is not a constant, taking a snapshot of its value" << std::endl;
92 }
93 (*_mref)[i] = static_cast<RooAbsReal *>(mref)->getVal();
94 ++i;
95 }
96
97 // initialization
98 initialize();
99}
100
101//_____________________________________________________________________________
104 _cacheMgr(other._cacheMgr, this),
105 m("m", this, other.m),
106 _varList("varList", this, other._varList),
107 _pdfList("pdfList", this, other._pdfList),
108 _mref(new TVectorD(*other._mref)),
109 _setting(other._setting),
110 _useHorizMorph(other._useHorizMorph)
111{
112
113 // initialization
114 initialize();
115}
116
117//_____________________________________________________________________________
119{
120 if (_mref)
121 delete _mref;
122 if (_M)
123 delete _M;
124}
125
126//_____________________________________________________________________________
128{
129
131
132 // other quantities needed
133 if (nPdf != _mref->GetNrows()) {
134 coutE(InputArguments) << "RooMomentMorphFunc::initialize(" << GetName() << ") ERROR: nPdf != nRefPoints" << std::endl;
135 assert(0);
136 }
137
138 TVectorD *dm = new TVectorD(nPdf);
139 _M = new TMatrixD(nPdf, nPdf);
140
141 // transformation matrix for non-linear extrapolation, needed in evaluate()
142 TMatrixD M(nPdf, nPdf);
143 for (Int_t i = 0; i < _mref->GetNrows(); ++i) {
144 (*dm)[i] = (*_mref)[i] - (*_mref)[0];
145 M(i, 0) = 1.;
146 if (i > 0)
147 M(0, i) = 0.;
148 }
149 for (Int_t i = 1; i < _mref->GetNrows(); ++i) {
150 for (Int_t j = 1; j < _mref->GetNrows(); ++j) {
151 M(i, j) = std::pow((*dm)[i], (double)j);
152 }
153 }
154 (*_M) = M.Invert();
155
156 delete dm;
157}
158
159//_____________________________________________________________________________
161{
162 auto cache = static_cast<CacheElem *>(_cacheMgr.getObj(nullptr, static_cast<RooArgSet const*>(nullptr)));
163 if (cache) {
164 return cache;
165 }
168
169 RooAbsReal *null = nullptr;
178
180
182
183 // fraction parameters
184 RooArgList coefList("coefList");
185 RooArgList coefList2("coefList2");
186 for (Int_t i = 0; i < 2 * nPdf; ++i) {
187 std::string fracName = Form("frac_%d", i);
188
189 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), 1.);
190
191 fracl.add(*frac); // to be set later
192 if (i < nPdf) {
193 coefList.add(*static_cast<RooRealVar *>(fracl.at(i)));
194 } else {
195 coefList2.add(*static_cast<RooRealVar *>(fracl.at(i)));
196 }
197 ownedComps.add(*static_cast<RooRealVar *>(fracl.at(i)));
198 }
199
200 RooRealSumFunc *theSumFunc = nullptr;
201 std::string sumfuncName = Form("%s_sumfunc", GetName());
202
203 if (_useHorizMorph) {
204 // mean and sigma
206 for (Int_t i = 0; i < nPdf; ++i) {
207 for (Int_t j = 0; j < nVar; ++j) {
208
209 std::string meanName = Form("%s_mean_%d_%d", GetName(), i, j);
210 std::string sigmaName = Form("%s_sigma_%d_%d", GetName(), i, j);
211
212 RooAbsMoment *mom = nVar == 1 ? (static_cast<RooAbsPdf *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*varList.at(j)))
213 : (static_cast<RooAbsPdf *>(_pdfList.at(i)))->sigma(static_cast<RooRealVar &>(*varList.at(j)), varList);
214
215 mom->setLocalNoDirtyInhibit(true);
216 mom->mean()->setLocalNoDirtyInhibit(true);
217
218 sigmarv[ij(i, j)] = mom;
219 meanrv[ij(i, j)] = mom->mean();
220
221 ownedComps.add(*sigmarv[ij(i, j)]);
222 }
223 }
224
225 // slope and offset (to be set later, depend on m)
226 for (Int_t j = 0; j < nVar; ++j) {
227 RooArgList meanList("meanList");
228 RooArgList rmsList("rmsList");
229 for (Int_t i = 0; i < nPdf; ++i) {
230 meanList.add(*meanrv[ij(i, j)]);
231 rmsList.add(*sigmarv[ij(i, j)]);
232 }
233 std::string myrmsName = Form("%s_rms_%d", GetName(), j);
234 std::string myposName = Form("%s_pos_%d", GetName(), j);
235 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList2);
236 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
237 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
238 }
239
240 // construction of unit pdfs
242
243 for (Int_t i = 0; i < nPdf; ++i) {
244 auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
245 std::string pdfName = Form("pdf_%d", i);
246 RooCustomizer cust(pdf, pdfName.c_str());
247
248 for (Int_t j = 0; j < nVar; ++j) {
249 // slope and offset formulas
250 std::string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
251 std::string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
252 slope[ij(i, j)] = new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[ij(i, j)], *myrms[j]));
253 offs[ij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
254 RooArgList(*meanrv[ij(i, j)], *mypos[j], *slope[ij(i, j)]));
255 ownedComps.add(RooArgSet(*slope[ij(i, j)], *offs[ij(i, j)]));
256 // linear transformations, so pdf can be renormalized
257 auto& var = static_cast<RooRealVar&>(*_varList[j]);
258 std::string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
259 // transVar[ij(i,j)] = new
260 // RooFormulaVar(transVarName.c_str(),transVarName.c_str(),"@0*@1+@2",RooArgList(*var,*slope[ij(i,j)],*offs[ij(i,j)]));
261
262 transVar[ij(i, j)] =
263 new RooLinearVar(transVarName.c_str(), transVarName.c_str(), var, *slope[ij(i, j)], *offs[ij(i, j)]);
264
265 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
266 // This will prevent the likelihood optimizers from erroneously declaring terms constant
267 transVar[ij(i, j)]->addServer((RooAbsArg &)m.arg());
268
269 ownedComps.add(*transVar[ij(i, j)]);
270 cust.replaceArg(var, *transVar[ij(i, j)]);
271 }
272 transPdf[i] = static_cast<RooAbsPdf *>(cust.build());
273 transPdfList.add(*transPdf[i]);
274 ownedComps.add(*transPdf[i]);
275 }
276 // sum pdf
277 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), transPdfList, coefList);
278 } else {
279 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), _pdfList, coefList);
280 }
281
282 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
283 // This will prevent the likelihood optimizers from erroneously declaring terms constant
284 theSumFunc->addServer((RooAbsArg &)m.arg());
285 theSumFunc->addOwnedComponents(ownedComps);
286
287 // change tracker for fraction parameters
288 std::string trackerName = Form("%s_frac_tracker", GetName());
289 RooChangeTracker *tracker = new RooChangeTracker(trackerName.c_str(), trackerName.c_str(), m.arg(), true);
290
291 // Store it in the cache
292 cache = new CacheElem(*theSumFunc, *tracker, fracl);
293 _cacheMgr.setObj(nullptr, nullptr, cache, nullptr);
294
295 return cache;
296}
297
298//_____________________________________________________________________________
303
304//_____________________________________________________________________________
306{
307 delete _sumFunc;
308 delete _tracker;
309}
310
311//_____________________________________________________________________________
313{
314 // Special version of getValV() overrides RooAbsReal::getVal() to save value of current normalization set
315 _curNormSet = set ? const_cast<RooArgSet *>(set) : const_cast<RooArgSet *>(static_cast<RooArgSet const*>(&_varList));
316 return RooAbsReal::getValV(set);
317}
318
319//_____________________________________________________________________________
321{
322 CacheElem *cache = getCache(nset ? nset : _curNormSet);
323
324 if (cache->_tracker->hasChanged(true)) {
325 cache->calculateFractions(*this, false); // verbose turned off
326 }
327
328 return cache->_sumFunc;
329}
330
331//_____________________________________________________________________________
333{
334 CacheElem *cache = getCache(nset ? nset : _curNormSet);
335
336 if (cache->_tracker->hasChanged(true)) {
337 cache->calculateFractions(*this, false); // verbose turned off
338 }
339
340 return cache->_sumFunc;
341}
342
343//_____________________________________________________________________________
345{
347
348 if (cache->_tracker->hasChanged(true)) {
349 cache->calculateFractions(*this, false); // verbose turned off
350 }
351
352 double ret = cache->_sumFunc->getVal(_pdfList.nset());
353 return ret;
354}
355
356//_____________________________________________________________________________
358{
359 return static_cast<RooRealVar *>(_frac.at(i));
360}
361
362//_____________________________________________________________________________
364{
365 return static_cast<RooRealVar *>(_frac.at(i));
366}
367
368//_____________________________________________________________________________
370{
371 Int_t nPdf = self._pdfList.size();
372
373 double dm = self.m - (*self._mref)[0];
374
375 // fully non-linear
376 double sumposfrac = 0.;
377 for (Int_t i = 0; i < nPdf; ++i) {
378 double ffrac = 0.;
379 for (Int_t j = 0; j < nPdf; ++j) {
380 ffrac += (*self._M)(j, i) * (j == 0 ? 1. : std::pow(dm, (double)j));
381 }
382 if (ffrac >= 0)
383 sumposfrac += ffrac;
384 // fractions for pdf
385 const_cast<RooRealVar *>(frac(i))->setVal(ffrac);
386 // fractions for rms and mean
387 const_cast<RooRealVar *>(frac(nPdf + i))->setVal(ffrac);
388 if (verbose) {
389 std::cout << ffrac << std::endl;
390 }
391 }
392
393 // various mode settings
394 int imin = self.idxmin(self.m);
395 int imax = self.idxmax(self.m);
396 double mfrac = (self.m - (*self._mref)[imin]) / ((*self._mref)[imax] - (*self._mref)[imin]);
397 switch (self._setting) {
398 case NonLinear:
399 // default already set above
400 break;
401
402 case SineLinear:
403 mfrac =
404 std::sin(TMath::PiOver2() * mfrac); // this gives a continuous differentiable transition between grid points.
405
406 // now fall through to Linear case
407
408 case Linear:
409 for (Int_t i = 0; i < 2 * nPdf; ++i) const_cast<RooRealVar *>(frac(i))->setVal(0.);
410 if (imax > imin) { // m in between mmin and mmax
411 const_cast<RooRealVar *>(frac(imin))->setVal(1. - mfrac);
412 const_cast<RooRealVar *>(frac(nPdf + imin))->setVal(1. - mfrac);
413 const_cast<RooRealVar *>(frac(imax))->setVal(mfrac);
414 const_cast<RooRealVar *>(frac(nPdf + imax))->setVal(mfrac);
415 } else if (imax == imin) { // m outside mmin and mmax
416 const_cast<RooRealVar *>(frac(imin))->setVal(1.);
417 const_cast<RooRealVar *>(frac(nPdf + imin))->setVal(1.);
418 }
419 break;
421 for (Int_t i = 0; i < nPdf; ++i) const_cast<RooRealVar *>(frac(i))->setVal(0.);
422 if (imax > imin) { // m in between mmin and mmax
423 const_cast<RooRealVar *>(frac(imin))->setVal(1. - mfrac);
424 const_cast<RooRealVar *>(frac(imax))->setVal(mfrac);
425 } else if (imax == imin) { // m outside mmin and mmax
426 const_cast<RooRealVar *>(frac(imin))->setVal(1.);
427 }
428 break;
430 for (Int_t i = 0; i < nPdf; ++i) {
431 if (frac(i)->getVal() < 0)
432 const_cast<RooRealVar *>(frac(i))->setVal(0.);
433 const_cast<RooRealVar *>(frac(i))->setVal(frac(i)->getVal() / sumposfrac);
434 }
435 break;
436 }
437}
438
439//_____________________________________________________________________________
440int RooMomentMorphFunc::idxmin(const double &mval) const
441{
442 int imin(0);
444 double mmin = -DBL_MAX;
445 for (Int_t i = 0; i < nPdf; ++i) {
446 if ((*_mref)[i] > mmin && (*_mref)[i] <= mval) {
447 mmin = (*_mref)[i];
448 imin = i;
449 }
450 }
451 return imin;
452}
453
454//_____________________________________________________________________________
455int RooMomentMorphFunc::idxmax(const double &mval) const
456{
457 int imax(0);
459 double mmax = DBL_MAX;
460 for (Int_t i = 0; i < nPdf; ++i) {
461 if ((*_mref)[i] < mmax && (*_mref)[i] >= mval) {
462 mmax = (*_mref)[i];
463 imax = i;
464 }
465 }
466 return imax;
467}
468
469//_____________________________________________________________________________
470std::list<double> *RooMomentMorphFunc::plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
471{
472 return sumFunc(nullptr)->plotSamplingHint(obs, xlo, xhi);
473}
474
475//_____________________________________________________________________________
476std::list<double> *RooMomentMorphFunc::binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
477{
478 return sumFunc(nullptr)->binBoundaries(obs, xlo, xhi);
479}
480
481//_____________________________________________________________________________
483{
484 return sumFunc(nullptr)->isBinnedDistribution(obs);
485}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutW(a)
#define coutE(a)
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
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
virtual std::list< double > * binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
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:402
virtual bool isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:348
virtual std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
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
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.
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.
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 RooMomentMorphFunc &self, bool verbose=true) const
RooArgList containedArgs(Action) override
double getValV(const RooArgSet *set=nullptr) const override
Return value of object.
friend class CacheElem
Current normalization set.
RooArgSet * _curNormSet
The cache manager.
bool isBinnedDistribution(const RooArgSet &obs) const override
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
RooObjCacheManager _cacheMgr
CacheElem * getCache(const RooArgSet *nset) const
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsReal * sumFunc(const RooArgSet *nset)
int idxmin(const double &m) const
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Retrieve bin boundaries if this distribution is binned in obs.
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
Int_t ij(const Int_t &i, const Int_t &j) const
int idxmax(const double &m) const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const T & arg() const
Return reference to object held in proxy.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Int_t GetNrows() const
Definition TVectorT.h:75
constexpr Double_t PiOver2()
Definition TMath.h:51
TMarker m
Definition textangle.C:8