Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorph.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 RooMomentMorph
12 \ingroup Roofit
13
14**/
15
16#include "Riostream.h"
17
18#include "RooMomentMorph.h"
19#include "RooAbsCategory.h"
20#include "RooRealConstant.h"
21#include "RooRealVar.h"
22#include "RooFormulaVar.h"
23#include "RooCustomizer.h"
24#include "RooAddPdf.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
33
34////////////////////////////////////////////////////////////////////////////////
35/// coverity[UNINIT_CTOR]
36
38 : _cacheMgr(this,10,true,true)
39{
40}
41
42////////////////////////////////////////////////////////////////////////////////
43/// CTOR
44
45RooMomentMorph::RooMomentMorph(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
46 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
47 : RooAbsPdf(name, title),
48 _cacheMgr(this, 10, true, true),
49 m("m", "m", this, _m),
50 _varList("varList", "List of variables", this),
51 _pdfList("pdfList", "List of pdfs", this),
52 _mref(new TVectorD(mrefpoints)),
53 _setting(setting),
54 _useHorizMorph(true)
55{
56 // observables
58
59 // reference p.d.f.s
60 _pdfList.addTyped<RooAbsPdf>(pdfList);
61
62 // initialization
63 initialize();
64}
65
66////////////////////////////////////////////////////////////////////////////////
67/// CTOR
68
69RooMomentMorph::RooMomentMorph(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
70 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
71 : RooAbsPdf(name, title),
72 _cacheMgr(this, 10, true, true),
73 m("m", "m", this, _m),
74 _varList("varList", "List of variables", this),
75 _pdfList("pdfList", "List of pdfs", this),
76 _mref(new TVectorD(mrefList.size())),
77 _setting(setting),
78 _useHorizMorph(true)
79{
80 // observables
82
83 // reference p.d.f.s
84 _pdfList.addTyped<RooAbsPdf>(pdfList);
85
86 // reference points in m
87
88 Int_t i = 0;
89 for (auto *mref : mrefList) {
90 if (!dynamic_cast<RooAbsReal*>(mref)) {
91 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: mref " << mref->GetName() << " is not of type RooAbsReal" << std::endl ;
92 throw std::string("RooPolyMorh::ctor() ERROR mref is not of type RooAbsReal") ;
93 }
94 if (!dynamic_cast<RooConstVar*>(mref)) {
95 coutW(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") WARNING mref point " << i << " is not a constant, taking a snapshot of its value" << std::endl ;
96 }
97 (*_mref)[i] = static_cast<RooAbsReal*>(mref)->getVal() ;
98 i++;
99 }
100
101 // initialization
102 initialize();
103}
104
105////////////////////////////////////////////////////////////////////////////////
106
108 : RooAbsPdf(other, name),
109 _cacheMgr(other._cacheMgr, this),
110 m("m", this, other.m),
111 _varList("varList", this, other._varList),
112 _pdfList("pdfList", this, other._pdfList),
113 _mref(new TVectorD(*other._mref)),
114 _setting(other._setting),
115 _useHorizMorph(other._useHorizMorph)
116{
117
118 // initialization
119 initialize();
120}
121
122////////////////////////////////////////////////////////////////////////////////
123
125{
126 if (_mref) delete _mref;
127 if (_M) delete _M;
128}
129
130////////////////////////////////////////////////////////////////////////////////
131
133{
135
136 // other quantities needed
137 if (nPdf!=_mref->GetNrows()) {
138 coutE(InputArguments) << "RooMomentMorph::initialize(" << GetName() << ") ERROR: nPdf != nRefPoints" << std::endl ;
139 assert(0) ;
140 }
141
142 TVectorD dm{nPdf};
143 _M = new TMatrixD(nPdf,nPdf);
144
145 // transformation matrix for non-linear extrapolation, needed in evaluate()
146 TMatrixD M(nPdf,nPdf);
147 for (Int_t i=0; i<_mref->GetNrows(); ++i) {
148 dm[i] = (*_mref)[i]-(*_mref)[0];
149 M(i,0) = 1.;
150 if (i>0) M(0,i) = 0.;
151 }
152 for (Int_t i=1; i<_mref->GetNrows(); ++i) {
153 for (Int_t j=1; j<_mref->GetNrows(); ++j) {
154 M(i,j) = std::pow(dm[i],(double)j);
155 }
156 }
157 (*_M) = M.Invert();
158}
159
160////////////////////////////////////////////////////////////////////////////////
161
162RooMomentMorph::CacheElem::CacheElem(std::unique_ptr<RooAbsPdf> && sumPdf,
163 std::unique_ptr<RooChangeTracker> && tracker,
164 const RooArgList& flist)
165 : _sumPdf(std::move(sumPdf)), _tracker(std::move(tracker)) {
166 _frac.add(flist);
167}
168
169////////////////////////////////////////////////////////////////////////////////
170
172{
173 if (auto* cache = static_cast<CacheElem*>(_cacheMgr.getObj(nullptr,static_cast<RooArgSet*>(nullptr)))) {
174 return cache ;
175 }
178
179 RooAbsReal* null = nullptr;
180 std::vector<RooAbsReal*> meanrv(nPdf*nVar,null);
181 std::vector<RooAbsReal*> sigmarv(nPdf*nVar,null);
182 std::vector<RooAbsReal*> myrms(nVar,null);
183 std::vector<RooAbsReal*> mypos(nVar,null);
184 std::vector<RooAbsReal*> slope(nPdf*nVar,null);
185 std::vector<RooAbsReal*> offs(nPdf*nVar,null);
186 std::vector<RooAbsReal*> transVar(nPdf*nVar,null);
187 std::vector<RooAbsReal*> transPdf(nPdf,null);
188
190
192
193 // fraction parameters
194 RooArgList coefList("coefList");
195 RooArgList coefList2("coefList2");
196 for (Int_t i=0; i<2*nPdf; ++i) {
197 std::string fracName = "frac_" + std::to_string(i);
198
199 RooRealVar* frac = new RooRealVar(fracName.c_str(),fracName.c_str(),1.) ;
200
201 fracl.add(*frac); // to be set later
202 if (i<nPdf) coefList.add(*frac);
203 else coefList2.add(*frac);
204 ownedComps.add(*frac);
205 }
206
207 std::unique_ptr<RooAddPdf> theSumPdf;
208 std::string sumpdfName = Form("%s_sumpdf",GetName());
209
210 if (_useHorizMorph) {
211 // mean and sigma
213 for (Int_t i=0; i<nPdf; ++i) {
214 for (Int_t j=0; j<nVar; ++j) {
215
216 std::string meanName = Form("%s_mean_%d_%d",GetName(),i,j);
217 std::string sigmaName = Form("%s_sigma_%d_%d",GetName(),i,j);
218
219 RooAbsMoment* mom = nVar==1 ?
220 (static_cast<RooAbsPdf*>(_pdfList.at(i)))->sigma(static_cast<RooRealVar&>(*varList.at(j))) :
221 (static_cast<RooAbsPdf*>(_pdfList.at(i)))->sigma(static_cast<RooRealVar&>(*varList.at(j)),varList) ;
222
223 mom->setLocalNoDirtyInhibit(true) ;
224 mom->mean()->setLocalNoDirtyInhibit(true) ;
225
226 sigmarv[ij(i,j)] = mom ;
227 meanrv[ij(i,j)] = mom->mean() ;
228
229 ownedComps.add(*sigmarv[ij(i,j)]) ;
230 }
231 }
232
233 // slope and offset (to be set later, depend on m)
234 for (Int_t j=0; j<nVar; ++j) {
235 RooArgList meanList("meanList");
236 RooArgList rmsList("rmsList");
237 for (Int_t i=0; i<nPdf; ++i) {
238 meanList.add(*meanrv[ij(i,j)]);
239 rmsList.add(*sigmarv[ij(i,j)]);
240 }
241 std::string myrmsName = Form("%s_rms_%d",GetName(),j);
242 std::string myposName = Form("%s_pos_%d",GetName(),j);
243 myrms[j] = new RooAddition(myrmsName.c_str(),myrmsName.c_str(),rmsList,coefList2);
244 mypos[j] = new RooAddition(myposName.c_str(),myposName.c_str(),meanList,coefList2);
245 ownedComps.add(RooArgSet(*myrms[j],*mypos[j])) ;
246 }
247
248 // construction of unit pdfs
250
251 for (Int_t i=0; i<nPdf; ++i) {
252
253 auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
254 std::string pdfName = "pdf_" + std::to_string(i);
255 RooCustomizer cust(pdf,pdfName.c_str());
256
257 for (Int_t j=0; j<nVar; ++j) {
258 // slope and offset formulas
259 std::string slopeName = Form("%s_slope_%d_%d",GetName(),i,j);
260 std::string offsetName = Form("%s_offset_%d_%d",GetName(),i,j);
261 slope[ij(i,j)] = new RooFormulaVar(slopeName.c_str(),"@0/@1", {*sigmarv[ij(i,j)],*myrms[j]});
262 offs[ij(i,j)] = new RooFormulaVar(offsetName.c_str(),"@0-(@1*@2)", {*meanrv[ij(i,j)],*mypos[j],*slope[ij(i,j)]});
263 ownedComps.add(RooArgSet(*slope[ij(i,j)],*offs[ij(i,j)])) ;
264 // linear transformations, so pdf can be renormalized
265 auto* var = static_cast<RooRealVar*>(_varList[j]);
266 std::string transVarName = Form("%s_transVar_%d_%d",GetName(),i,j);
267 //transVar[ij(i,j)] = new RooFormulaVar(transVarName.c_str(),transVarName.c_str(),"@0*@1+@2",RooArgList(*var,*slope[ij(i,j)],*offs[ij(i,j)]));
268
269 transVar[ij(i,j)] = new RooLinearVar(transVarName.c_str(),transVarName.c_str(),*var,*slope[ij(i,j)],*offs[ij(i,j)]);
270
271 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
272 // This will prevent the likelihood optimizers from erroneously declaring terms constant
273 transVar[ij(i,j)]->addServer(const_cast<RooAbsReal&>(m.arg()));
274
275 ownedComps.add(*transVar[ij(i,j)]) ;
276 cust.replaceArg(*var,*transVar[ij(i,j)]);
277 }
278 transPdf[i] = static_cast<RooAbsPdf*>(cust.build()) ;
279 transPdfList.add(*transPdf[i]);
280 ownedComps.add(*transPdf[i]) ;
281 }
282 // sum pdf
283 theSumPdf = std::make_unique<RooAddPdf>(sumpdfName.c_str(),sumpdfName.c_str(),transPdfList,coefList);
284 }
285 else {
286 theSumPdf = std::make_unique<RooAddPdf>(sumpdfName.c_str(),sumpdfName.c_str(),_pdfList,coefList);
287 }
288 theSumPdf->fixCoefNormalization(*nset);
289
290 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
291 // This will prevent the likelihood optimizers from erroneously declaring terms constant
292 theSumPdf->addServer(const_cast<RooAbsReal&>(m.arg()));
293 theSumPdf->addOwnedComponents(ownedComps) ;
294
295 // change tracker for fraction parameters
296 std::string trackerName = Form("%s_frac_tracker",GetName()) ;
297 auto tracker = std::make_unique<RooChangeTracker>(trackerName.c_str(),trackerName.c_str(),m.arg(),true) ;
298
299 // Store it in the cache
300 auto cache = new CacheElem(std::move(theSumPdf),std::move(tracker),fracl) ;
301 _cacheMgr.setObj(nullptr,nullptr,cache,nullptr);
302
303 cache->calculateFractions(*this, false);
304 return cache ;
305}
306
307////////////////////////////////////////////////////////////////////////////////
308
310{
311 return RooArgList(*_sumPdf,*_tracker) ;
312}
313
314////////////////////////////////////////////////////////////////////////////////
315
317
318////////////////////////////////////////////////////////////////////////////////
319/// Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
320
321double RooMomentMorph::getValV(const RooArgSet* set) const
322{
323 _curNormSet = set ? const_cast<RooArgSet*>(set) : const_cast<RooArgSet*>(static_cast<RooArgSet const*>(&_varList));
324 return RooAbsPdf::getValV(set) ;
325}
326
327////////////////////////////////////////////////////////////////////////////////
328
330{
331 CacheElem* cache = getCache(nset ? nset : _curNormSet) ;
332
333 if (cache->_tracker->hasChanged(true)) {
334 cache->calculateFractions(*this,false); // verbose turned off
335 }
336
337 return cache->_sumPdf.get();
338}
339
340////////////////////////////////////////////////////////////////////////////////
341
343{
344 CacheElem* cache = getCache(_curNormSet) ;
345
346 if (cache->_tracker->hasChanged(true)) {
347 cache->calculateFractions(*this,false); // verbose turned off
348 }
349
350 return cache->_sumPdf->getVal(_pdfList.nset());
351}
352
353////////////////////////////////////////////////////////////////////////////////
354
356{
357 return static_cast<RooRealVar*>(_frac.at(i)) ;
358}
359
360////////////////////////////////////////////////////////////////////////////////
361
363{
364 return static_cast<RooRealVar*>(_frac.at(i)) ;
365}
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) { ffrac += (*self._M)(j,i) * (j==0?1.:std::pow(dm,(double)j)); }
380 if (ffrac>=0) sumposfrac+=ffrac;
381 // fractions for pdf
382 const_cast<RooRealVar*>(frac(i))->setVal(ffrac);
383 // fractions for rms and mean
384 const_cast<RooRealVar*>(frac(nPdf+i))->setVal(ffrac);
385 if (verbose) { std::cout << ffrac << std::endl; }
386 }
387
388 // various mode settings
389 int imin = self.idxmin(self.m);
390 int imax = self.idxmax(self.m);
391 double mfrac = (self.m-(*self._mref)[imin])/((*self._mref)[imax]-(*self._mref)[imin]);
392 switch (self._setting) {
393 case NonLinear:
394 // default already set above
395 break;
396
397 case SineLinear:
398 mfrac = std::sin( TMath::PiOver2()*mfrac ); // this gives a continuous differentiable transition between grid points.
399
400 // now fall through to Linear case
401
402 case Linear:
403 for (Int_t i=0; i<2*nPdf; ++i)
404 const_cast<RooRealVar*>(frac(i))->setVal(0.);
405 if (imax>imin) { // m in between mmin and mmax
406 const_cast<RooRealVar*>(frac(imin))->setVal(1.-mfrac);
407 const_cast<RooRealVar*>(frac(nPdf+imin))->setVal(1.-mfrac);
408 const_cast<RooRealVar*>(frac(imax))->setVal(mfrac);
409 const_cast<RooRealVar*>(frac(nPdf+imax))->setVal(mfrac);
410 } else if (imax==imin) { // m outside mmin and mmax
411 const_cast<RooRealVar*>(frac(imin))->setVal(1.);
412 const_cast<RooRealVar*>(frac(nPdf+imin))->setVal(1.);
413 }
414 break;
416 for (Int_t i=0; i<nPdf; ++i)
417 const_cast<RooRealVar*>(frac(i))->setVal(0.);
418 if (imax>imin) { // m in between mmin and mmax
419 const_cast<RooRealVar*>(frac(imin))->setVal(1.-mfrac);
420 const_cast<RooRealVar*>(frac(imax))->setVal(mfrac);
421 } else if (imax==imin) { // m outside mmin and mmax
422 const_cast<RooRealVar*>(frac(imin))->setVal(1.);
423 }
424 break;
426 for (Int_t i = 0; i < nPdf; ++i) {
427 if (frac(i)->getVal() < 0)
428 const_cast<RooRealVar *>(frac(i))->setVal(0.);
429 const_cast<RooRealVar *>(frac(i))->setVal(frac(i)->getVal() / sumposfrac);
430 }
431 break;
432 }
433
434}
435
436////////////////////////////////////////////////////////////////////////////////
437
438int RooMomentMorph::idxmin(const double& mval) const
439{
440 int imin(0);
442 double mmin=-DBL_MAX;
443 for (Int_t i=0; i<nPdf; ++i)
444 if ( (*_mref)[i]>mmin && (*_mref)[i]<=mval ) { mmin=(*_mref)[i]; imin=i; }
445 return imin;
446}
447
448
449////////////////////////////////////////////////////////////////////////////////
450
451int RooMomentMorph::idxmax(const double& mval) const
452{
453 int imax(0);
455 double mmax=DBL_MAX;
456 for (Int_t i=0; i<nPdf; ++i)
457 if ( (*_mref)[i]<mmax && (*_mref)[i]>=mval ) { mmax=(*_mref)[i]; imax=i; }
458 return imax;
459}
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
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
double getValV(const RooArgSet *set=nullptr) const override
Return current value, normalized by integrating over the observables in nset.
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
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.
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...
RooArgList containedArgs(Action) override
CacheElem(std::unique_ptr< RooAbsPdf > &&sumPdf, std::unique_ptr< RooChangeTracker > &&tracker, const RooArgList &flist)
RooRealVar * frac(Int_t i)
void calculateFractions(const RooMomentMorph &self, bool verbose=true) const
std::unique_ptr< RooAbsPdf > _sumPdf
std::unique_ptr< RooChangeTracker > _tracker
RooObjCacheManager _cacheMgr
friend class CacheElem
Current normalization set.
~RooMomentMorph() override
RooSetProxy _varList
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooMomentMorph()
coverity[UNINIT_CTOR]
RooRealProxy m
RooAbsPdf * sumPdf(const RooArgSet *nset)
Int_t ij(const Int_t &i, const Int_t &j) const
RooArgSet * _curNormSet
The cache manager.
RooListProxy _pdfList
int idxmin(const double &m) const
int idxmax(const double &m) const
double getValV(const RooArgSet *set=nullptr) const override
Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set...
CacheElem * getCache(const RooArgSet *nset) const
TVectorD * _mref
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
const T & arg() const
Return reference to object held in proxy.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Int_t GetNrows() const
Definition TVectorT.h:75
constexpr Double_t PiOver2()
Definition TMath.h:51
TMarker m
Definition textangle.C:8