Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorph.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * This code was autogenerated by RooClassFactory *
5 *****************************************************************************/
6
7/** \class RooMomentMorph
8 \ingroup Roofit
9
10**/
11
12#include "Riostream.h"
13
14#include "RooMomentMorph.h"
15#include "RooAbsCategory.h"
16#include "RooRealIntegral.h"
17#include "RooRealConstant.h"
18#include "RooRealVar.h"
19#include "RooFormulaVar.h"
20#include "RooCustomizer.h"
21#include "RooAddPdf.h"
22#include "RooAddition.h"
23#include "RooMoment.h"
24#include "RooLinearVar.h"
25#include "RooChangeTracker.h"
26
27#include "TMath.h"
28#include "TH1.h"
29
30using namespace std;
31
33
34////////////////////////////////////////////////////////////////////////////////
35/// coverity[UNINIT_CTOR]
36
38 : _cacheMgr(this,10,true,true), _curNormSet(0), _mref(0), _M(0), _useHorizMorph(true)
39{
42}
43
44////////////////////////////////////////////////////////////////////////////////
45/// CTOR
46
47RooMomentMorph::RooMomentMorph(const char *name, const char *title,
48 RooAbsReal& _m,
49 const RooArgList& varList,
50 const RooArgList& pdfList,
51 const TVectorD& mrefpoints,
52 Setting setting) :
53 RooAbsPdf(name,title),
54 _cacheMgr(this,10,kTRUE,kTRUE),
55 m("m","m",this,_m),
56 _varList("varList","List of variables",this),
57 _pdfList("pdfList","List of pdfs",this),
58 _setting(setting),
59 _useHorizMorph(true)
60{
61 // observables
62 TIterator* varItr = varList.createIterator() ;
63 RooAbsArg* var ;
64 for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
65 if (!dynamic_cast<RooAbsReal*>(var)) {
66 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: variable " << var->GetName() << " is not of type RooAbsReal" << endl ;
67 throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal") ;
68 }
69 _varList.add(*var) ;
70 }
71 delete varItr ;
72
73 // reference p.d.f.s
74 TIterator* pdfItr = pdfList.createIterator() ;
75 RooAbsPdf* pdf ;
76 for (Int_t i=0; (pdf = dynamic_cast<RooAbsPdf*>(pdfItr->Next())); ++i) {
77 if (!pdf) {
78 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: pdf " << pdf->GetName() << " is not of type RooAbsPdf" << endl ;
79 throw string("RooPolyMorh::ctor() ERROR pdf is not of type RooAbsPdf") ;
80 }
81 _pdfList.add(*pdf) ;
82 }
83 delete pdfItr ;
84
85 _mref = new TVectorD(mrefpoints);
88
89 // initialization
90 initialize();
91}
92
93////////////////////////////////////////////////////////////////////////////////
94/// CTOR
95
96RooMomentMorph::RooMomentMorph(const char *name, const char *title,
97 RooAbsReal& _m,
98 const RooArgList& varList,
99 const RooArgList& pdfList,
100 const RooArgList& mrefList,
101 Setting setting) :
102 RooAbsPdf(name,title),
103 _cacheMgr(this,10,kTRUE,kTRUE),
104 m("m","m",this,_m),
105 _varList("varList","List of variables",this),
106 _pdfList("pdfList","List of pdfs",this),
107 _setting(setting),
108 _useHorizMorph(true)
109{
110 // observables
111 TIterator* varItr = varList.createIterator() ;
112 RooAbsArg* var ;
113 for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
114 if (!dynamic_cast<RooAbsReal*>(var)) {
115 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: variable " << var->GetName() << " is not of type RooAbsReal" << endl ;
116 throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal") ;
117 }
118 _varList.add(*var) ;
119 }
120 delete varItr ;
121
122 // reference p.d.f.s
123 TIterator* pdfItr = pdfList.createIterator() ;
124 RooAbsPdf* pdf ;
125 for (Int_t i=0; (pdf = dynamic_cast<RooAbsPdf*>(pdfItr->Next())); ++i) {
126 if (!pdf) {
127 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: pdf " << pdf->GetName() << " is not of type RooAbsPdf" << endl ;
128 throw string("RooPolyMorh::ctor() ERROR pdf is not of type RooAbsPdf") ;
129 }
130 _pdfList.add(*pdf) ;
131 }
132 delete pdfItr ;
133
134 // reference points in m
135 _mref = new TVectorD(mrefList.getSize());
136 TIterator* mrefItr = mrefList.createIterator() ;
137 RooAbsReal* mref ;
138 for (Int_t i=0; (mref = dynamic_cast<RooAbsReal*>(mrefItr->Next())); ++i) {
139 if (!mref) {
140 coutE(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") ERROR: mref " << mref->GetName() << " is not of type RooAbsReal" << endl ;
141 throw string("RooPolyMorh::ctor() ERROR mref is not of type RooAbsReal") ;
142 }
143 if (!dynamic_cast<RooConstVar*>(mref)) {
144 coutW(InputArguments) << "RooMomentMorph::ctor(" << GetName() << ") WARNING mref point " << i << " is not a constant, taking a snapshot of its value" << endl ;
145 }
146 (*_mref)[i] = mref->getVal() ;
147 }
148 delete mrefItr ;
149
152
153 // initialization
154 initialize();
155}
156
157////////////////////////////////////////////////////////////////////////////////
158
160 RooAbsPdf(other,name),
161 _cacheMgr(other._cacheMgr,this),
162 _curNormSet(0),
163 m("m",this,other.m),
164 _varList("varList",this,other._varList),
165 _pdfList("pdfList",this,other._pdfList),
166 _setting(other._setting),
167 _useHorizMorph(other._useHorizMorph)
168{
169 _mref = new TVectorD(*other._mref) ;
172
173 // initialization
174 initialize();
175}
176
177////////////////////////////////////////////////////////////////////////////////
178
180{
181 if (_mref) delete _mref;
182 if (_varItr) delete _varItr;
183 if (_pdfItr) delete _pdfItr;
184 if (_M) delete _M;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188
190{
191 Int_t nPdf = _pdfList.getSize();
192
193 // other quantities needed
194 if (nPdf!=_mref->GetNrows()) {
195 coutE(InputArguments) << "RooMomentMorph::initialize(" << GetName() << ") ERROR: nPdf != nRefPoints" << endl ;
196 assert(0) ;
197 }
198
199 TVectorD* dm = new TVectorD(nPdf);
200 _M = new TMatrixD(nPdf,nPdf);
201
202 // transformation matrix for non-linear extrapolation, needed in evaluate()
203 TMatrixD M(nPdf,nPdf);
204 for (Int_t i=0; i<_mref->GetNrows(); ++i) {
205 (*dm)[i] = (*_mref)[i]-(*_mref)[0];
206 M(i,0) = 1.;
207 if (i>0) M(0,i) = 0.;
208 }
209 for (Int_t i=1; i<_mref->GetNrows(); ++i) {
210 for (Int_t j=1; j<_mref->GetNrows(); ++j) {
211 M(i,j) = TMath::Power((*dm)[i],(double)j);
212 }
213 }
214 (*_M) = M.Invert();
215
216 delete dm ;
217}
218
219////////////////////////////////////////////////////////////////////////////////
220
222{
223 CacheElem* cache = (CacheElem*) _cacheMgr.getObj(0,(RooArgSet*)0) ;
224 if (cache) {
225 return cache ;
226 }
227 Int_t nVar = _varList.getSize();
228 Int_t nPdf = _pdfList.getSize();
229
230 RooAbsReal* null = 0 ;
231 vector<RooAbsReal*> meanrv(nPdf*nVar,null);
232 vector<RooAbsReal*> sigmarv(nPdf*nVar,null);
233 vector<RooAbsReal*> myrms(nVar,null);
234 vector<RooAbsReal*> mypos(nVar,null);
235 vector<RooAbsReal*> slope(nPdf*nVar,null);
236 vector<RooAbsReal*> offs(nPdf*nVar,null);
237 vector<RooAbsReal*> transVar(nPdf*nVar,null);
238 vector<RooAbsReal*> transPdf(nPdf,null);
239
240 RooArgSet ownedComps ;
241
242 RooArgList fracl ;
243
244 // fraction parameters
245 RooArgList coefList("coefList");
246 RooArgList coefList2("coefList2");
247 for (Int_t i=0; i<2*nPdf; ++i) {
248 std::string fracName = Form("frac_%d",i);
249
250 RooRealVar* frac = new RooRealVar(fracName.c_str(),fracName.c_str(),1.) ;
251
252 fracl.add(*frac); // to be set later
253 if (i<nPdf) coefList.add(*(RooRealVar*)(fracl.at(i))) ;
254 else coefList2.add(*(RooRealVar*)(fracl.at(i))) ;
255 ownedComps.add(*(RooRealVar*)(fracl.at(i))) ;
256 }
257
258 RooAddPdf* theSumPdf = 0;
259 std::string sumpdfName = Form("%s_sumpdf",GetName());
260
261 if (_useHorizMorph) {
262 // mean and sigma
263 RooArgList varList(_varList) ;
264 for (Int_t i=0; i<nPdf; ++i) {
265 for (Int_t j=0; j<nVar; ++j) {
266
267 std::string meanName = Form("%s_mean_%d_%d",GetName(),i,j);
268 std::string sigmaName = Form("%s_sigma_%d_%d",GetName(),i,j);
269
270 RooAbsMoment* mom = nVar==1 ?
271 ((RooAbsPdf*)_pdfList.at(i))->sigma((RooRealVar&)*varList.at(j)) :
272 ((RooAbsPdf*)_pdfList.at(i))->sigma((RooRealVar&)*varList.at(j),varList) ;
273
276
277 sigmarv[ij(i,j)] = mom ;
278 meanrv[ij(i,j)] = mom->mean() ;
279
280 ownedComps.add(*sigmarv[ij(i,j)]) ;
281 }
282 }
283
284 // slope and offset (to be set later, depend on m)
285 for (Int_t j=0; j<nVar; ++j) {
286 RooArgList meanList("meanList");
287 RooArgList rmsList("rmsList");
288 for (Int_t i=0; i<nPdf; ++i) {
289 meanList.add(*meanrv[ij(i,j)]);
290 rmsList.add(*sigmarv[ij(i,j)]);
291 }
292 std::string myrmsName = Form("%s_rms_%d",GetName(),j);
293 std::string myposName = Form("%s_pos_%d",GetName(),j);
294 myrms[j] = new RooAddition(myrmsName.c_str(),myrmsName.c_str(),rmsList,coefList2);
295 mypos[j] = new RooAddition(myposName.c_str(),myposName.c_str(),meanList,coefList2);
296 ownedComps.add(RooArgSet(*myrms[j],*mypos[j])) ;
297 }
298
299 // construction of unit pdfs
300 _pdfItr->Reset();
301 RooAbsPdf* pdf;
302 RooArgList transPdfList;
303
304 for (Int_t i=0; i<nPdf; ++i) {
305 _varItr->Reset() ;
306 RooRealVar* var ;
307
308 pdf = (RooAbsPdf*)_pdfItr->Next();
309 std::string pdfName = Form("pdf_%d",i);
310 RooCustomizer cust(*pdf,pdfName.c_str());
311
312 for (Int_t j=0; j<nVar; ++j) {
313 // slope and offset formulas
314 std::string slopeName = Form("%s_slope_%d_%d",GetName(),i,j);
315 std::string offsetName = Form("%s_offset_%d_%d",GetName(),i,j);
316 slope[ij(i,j)] = new RooFormulaVar(slopeName.c_str(),"@0/@1",RooArgList(*sigmarv[ij(i,j)],*myrms[j]));
317 offs[ij(i,j)] = new RooFormulaVar(offsetName.c_str(),"@0-(@1*@2)",RooArgList(*meanrv[ij(i,j)],*mypos[j],*slope[ij(i,j)]));
318 ownedComps.add(RooArgSet(*slope[ij(i,j)],*offs[ij(i,j)])) ;
319 // linear transformations, so pdf can be renormalized
320 var = (RooRealVar*)(_varItr->Next());
321 std::string transVarName = Form("%s_transVar_%d_%d",GetName(),i,j);
322 //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)]));
323
324 transVar[ij(i,j)] = new RooLinearVar(transVarName.c_str(),transVarName.c_str(),*var,*slope[ij(i,j)],*offs[ij(i,j)]);
325
326 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
327 // This will prevent the likelihood optimizers from erroneously declaring terms constant
328 transVar[ij(i,j)]->addServer((RooAbsArg&)m.arg()) ;
329
330 ownedComps.add(*transVar[ij(i,j)]) ;
331 cust.replaceArg(*var,*transVar[ij(i,j)]);
332 }
333 transPdf[i] = (RooAbsPdf*) cust.build() ;
334 transPdfList.add(*transPdf[i]);
335 ownedComps.add(*transPdf[i]) ;
336 }
337 // sum pdf
338 theSumPdf = new RooAddPdf(sumpdfName.c_str(),sumpdfName.c_str(),transPdfList,coefList);
339 }
340 else {
341 theSumPdf = new RooAddPdf(sumpdfName.c_str(),sumpdfName.c_str(),_pdfList,coefList);
342 }
343
344 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
345 // This will prevent the likelihood optimizers from erroneously declaring terms constant
346 theSumPdf->addServer((RooAbsArg&)m.arg()) ;
347 theSumPdf->addOwnedComponents(ownedComps) ;
348
349 // change tracker for fraction parameters
350 std::string trackerName = Form("%s_frac_tracker",GetName()) ;
351 RooChangeTracker* tracker = new RooChangeTracker(trackerName.c_str(),trackerName.c_str(),m.arg(),kTRUE) ;
352
353 // Store it in the cache
354 cache = new CacheElem(*theSumPdf,*tracker,fracl) ;
355 _cacheMgr.setObj(0,0,cache,0) ;
356
357 cache->calculateFractions(*this, kFALSE);
358 return cache ;
359}
360
361////////////////////////////////////////////////////////////////////////////////
362
364{
365 return RooArgList(*_sumPdf,*_tracker) ;
366}
367
368////////////////////////////////////////////////////////////////////////////////
369
371{
372 delete _sumPdf ;
373 delete _tracker ;
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
378
380{
381 _curNormSet = set ? (RooArgSet*)set : (RooArgSet*)&_varList ;
382 return RooAbsPdf::getVal(set) ;
383}
384
385////////////////////////////////////////////////////////////////////////////////
386
388{
389 CacheElem* cache = getCache(nset ? nset : _curNormSet) ;
390
391 if (cache->_tracker->hasChanged(kTRUE)) {
392 cache->calculateFractions(*this,kFALSE); // verbose turned off
393 }
394
395 return cache->_sumPdf ;
396}
397
398////////////////////////////////////////////////////////////////////////////////
399
401{
402 CacheElem* cache = getCache(_curNormSet) ;
403
404 if (cache->_tracker->hasChanged(kTRUE)) {
405 cache->calculateFractions(*this,kFALSE); // verbose turned off
406 }
407
408 Double_t ret = cache->_sumPdf->getVal(_pdfList.nset());
409 return ret ;
410}
411
412////////////////////////////////////////////////////////////////////////////////
413
415{
416 return (RooRealVar*)(_frac.at(i)) ;
417}
418
419////////////////////////////////////////////////////////////////////////////////
420
422{
423 return (RooRealVar*)(_frac.at(i)) ;
424}
425
426////////////////////////////////////////////////////////////////////////////////
427
429{
430 Int_t nPdf = self._pdfList.getSize();
431
432 Double_t dm = self.m - (*self._mref)[0];
433
434 // fully non-linear
435 double sumposfrac=0.;
436 for (Int_t i=0; i<nPdf; ++i) {
437 double ffrac=0.;
438 for (Int_t j=0; j<nPdf; ++j) { ffrac += (*self._M)(j,i) * (j==0?1.:TMath::Power(dm,(double)j)); }
439 if (ffrac>=0) sumposfrac+=ffrac;
440 // fractions for pdf
441 ((RooRealVar*)frac(i))->setVal(ffrac);
442 // fractions for rms and mean
443 ((RooRealVar*)frac(nPdf+i))->setVal(ffrac);
444 if (verbose) { cout << ffrac << endl; }
445 }
446
447 // various mode settings
448 int imin = self.idxmin(self.m);
449 int imax = self.idxmax(self.m);
450 double mfrac = (self.m-(*self._mref)[imin])/((*self._mref)[imax]-(*self._mref)[imin]);
451 switch (self._setting) {
452 case NonLinear:
453 // default already set above
454 break;
455
456 case SineLinear:
457 mfrac = TMath::Sin( TMath::PiOver2()*mfrac ); // this gives a continuous differentiable transition between grid points.
458
459 // now fall through to Linear case
460
461 case Linear:
462 for (Int_t i=0; i<2*nPdf; ++i)
463 ((RooRealVar*)frac(i))->setVal(0.);
464 if (imax>imin) { // m in between mmin and mmax
465 ((RooRealVar*)frac(imin))->setVal(1.-mfrac);
466 ((RooRealVar*)frac(nPdf+imin))->setVal(1.-mfrac);
467 ((RooRealVar*)frac(imax))->setVal(mfrac);
468 ((RooRealVar*)frac(nPdf+imax))->setVal(mfrac);
469 } else if (imax==imin) { // m outside mmin and mmax
470 ((RooRealVar*)frac(imin))->setVal(1.);
471 ((RooRealVar*)frac(nPdf+imin))->setVal(1.);
472 }
473 break;
475 for (Int_t i=0; i<nPdf; ++i)
476 ((RooRealVar*)frac(i))->setVal(0.);
477 if (imax>imin) { // m in between mmin and mmax
478 ((RooRealVar*)frac(imin))->setVal(1.-mfrac);
479 ((RooRealVar*)frac(imax))->setVal(mfrac);
480 } else if (imax==imin) { // m outside mmin and mmax
481 ((RooRealVar*)frac(imin))->setVal(1.);
482 }
483 break;
485 for (Int_t i=0; i<nPdf; ++i) {
486 if (((RooRealVar*)frac(i))->getVal()<0) ((RooRealVar*)frac(i))->setVal(0.);
487 ((RooRealVar*)frac(i))->setVal(((RooRealVar*)frac(i))->getVal()/sumposfrac);
488 }
489 break;
490 }
491
492}
493
494////////////////////////////////////////////////////////////////////////////////
495
496int RooMomentMorph::idxmin(const double& mval) const
497{
498 int imin(0);
499 Int_t nPdf = _pdfList.getSize();
500 double mmin=-DBL_MAX;
501 for (Int_t i=0; i<nPdf; ++i)
502 if ( (*_mref)[i]>mmin && (*_mref)[i]<=mval ) { mmin=(*_mref)[i]; imin=i; }
503 return imin;
504}
505
506
507////////////////////////////////////////////////////////////////////////////////
508
509int RooMomentMorph::idxmax(const double& mval) const
510{
511 int imax(0);
512 Int_t nPdf = _pdfList.getSize();
513 double mmax=DBL_MAX;
514 for (Int_t i=0; i<nPdf; ++i)
515 if ( (*_mref)[i]<mmax && (*_mref)[i]>=mval ) { mmax=(*_mref)[i]; imax=i; }
516 return imax;
517}
#define coutW(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
TMatrixT< Double_t > TMatrixD
Definition TMatrixDfwd.h:23
char * Form(const char *fmt,...)
TVectorT< Double_t > TVectorD
Definition TVectorDfwd.h:23
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
friend class RooArgSet
Definition RooAbsArg.h:642
void addServer(RooAbsArg &server, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
void setLocalNoDirtyInhibit(Bool_t flag) const
Definition RooAbsArg.h:720
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
RooAbsReal * mean()
const RooArgSet * nset() const
Definition RooAbsProxy.h:45
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
friend class RooAddPdf
Definition RooAbsReal.h:558
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition RooAddition.h:27
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
Bool_t hasChanged(Bool_t clearState)
Returns true if state has changed since last call with clearState=kTRUE.
RooConstVar represent a constant real-valued object.
Definition RooConstVar.h:26
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
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...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
void calculateFractions(const RooMomentMorph &self, Bool_t verbose=kTRUE) const
RooRealVar * frac(Int_t i)
virtual RooArgList containedArgs(Action)
RooChangeTracker * _tracker
RooObjCacheManager _cacheMgr
TIterator * _pdfItr
do not persist
RooSetProxy _varList
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
virtual ~RooMomentMorph()
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
TIterator * _varItr
int idxmin(const double &m) const
virtual Double_t getVal(const RooArgSet *set=0) const
Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set...
int idxmax(const double &m) const
CacheElem * getCache(const RooArgSet *nset) const
TVectorD * _mref
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
const T & arg() const
Return reference to object held in proxy.
Iterator abstract base class.
Definition TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Int_t GetNrows() const
Definition TVectorT.h:75
const Double_t sigma
constexpr Double_t PiOver2()
Definition TMath.h:51
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:685
Double_t Sin(Double_t)
Definition TMath.h:589
auto * m
Definition textangle.C:8