Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooMomentMorphFunc.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * This code was autogenerated by RooClassFactory *
5 *****************************************************************************/
6
7// Your description goes here...
8
9#include "Riostream.h"
10
11#include "RooMomentMorphFunc.h"
12#include "RooAbsCategory.h"
13#include "RooRealIntegral.h"
14#include "RooRealConstant.h"
15#include "RooRealVar.h"
16#include "RooFormulaVar.h"
17#include "RooCustomizer.h"
18#include "RooRealSumFunc.h"
19#include "RooAddition.h"
20#include "RooMoment.h"
21#include "RooLinearVar.h"
22#include "RooChangeTracker.h"
23
24#include "TMath.h"
25#include "TH1.h"
26
27using namespace std;
28
30
31 //_____________________________________________________________________________
33 : _curNormSet(0), _mref(0), _M(0), _useHorizMorph(true)
34{
35 // coverity[UNINIT_CTOR]
36 _varItr = _varList.createIterator();
37 _pdfItr = _pdfList.createIterator();
38}
39
40//_____________________________________________________________________________
41RooMomentMorphFunc::RooMomentMorphFunc(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
42 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
43 : RooAbsReal(name, title), _cacheMgr(this, 10, kTRUE, kTRUE), m("m", "m", this, _m),
44 _varList("varList", "List of variables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
45 _useHorizMorph(true)
46{
47 // CTOR
48
49 // observables
50 TIterator *varItr = varList.createIterator();
51 RooAbsArg *var;
52 for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
53 if (!dynamic_cast<RooAbsReal *>(var)) {
54 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
55 << " is not of type RooAbsReal" << endl;
56 throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal");
57 }
58 _varList.add(*var);
59 }
60 delete varItr;
61
62 // reference p.d.f.s
63 TIterator *pdfItr = pdfList.createIterator();
64 RooAbsReal *pdf;
65 for (Int_t i = 0; (pdf = dynamic_cast<RooAbsReal *>(pdfItr->Next())); ++i) {
66 if (!pdf) {
67 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: func " << pdf->GetName()
68 << " is not of type RooAbsReal" << endl;
69 throw string("RooMomentMorhFunc::ctor() ERROR func is not of type RooAbsReal");
70 }
71 _pdfList.add(*pdf);
72 }
73 delete pdfItr;
74
75 _mref = new TVectorD(mrefpoints);
78
79 // initialization
80 initialize();
81}
82
83//_____________________________________________________________________________
84RooMomentMorphFunc::RooMomentMorphFunc(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
85 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
86 : RooAbsReal(name, title), _cacheMgr(this, 10, kTRUE, kTRUE), m("m", "m", this, _m),
87 _varList("varList", "List of variables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
88 _useHorizMorph(true)
89{
90 // CTOR
91
92 // observables
93 TIterator *varItr = varList.createIterator();
94 RooAbsArg *var;
95 for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
96 if (!dynamic_cast<RooAbsReal *>(var)) {
97 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
98 << " is not of type RooAbsReal" << endl;
99 throw string("RooMomentMorh::ctor() ERROR variable is not of type RooAbsReal");
100 }
101 _varList.add(*var);
102 }
103 delete varItr;
104
105 // reference p.d.f.s
106 TIterator *pdfItr = pdfList.createIterator();
107 RooAbsReal *pdf;
108 for (Int_t i = 0; (pdf = dynamic_cast<RooAbsReal *>(pdfItr->Next())); ++i) {
109 if (!pdf) {
110 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: function " << pdf->GetName()
111 << " is not of type RooAbsReal" << endl;
112 throw string("RooMomentMorh::ctor() ERROR function is not of type RooAbsReal");
113 }
114 _pdfList.add(*pdf);
115 }
116 delete pdfItr;
117
118 // reference points in m
119 _mref = new TVectorD(mrefList.getSize());
120 TIterator *mrefItr = mrefList.createIterator();
121 RooAbsReal *mref;
122 for (Int_t i = 0; (mref = dynamic_cast<RooAbsReal *>(mrefItr->Next())); ++i) {
123 if (!mref) {
124 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
125 << " is not of type RooAbsReal" << endl;
126 throw string("RooPolyMorh::ctor() ERROR mref is not of type RooAbsReal");
127 }
128 if (!dynamic_cast<RooConstVar *>(mref)) {
129 coutW(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") WARNING mref point " << i
130 << " is not a constant, taking a snapshot of its value" << endl;
131 }
132 (*_mref)[i] = mref->getVal();
133 }
134 delete mrefItr;
135
138
139 // initialization
140 initialize();
141}
142
143//_____________________________________________________________________________
145 : RooAbsReal(other, name), _cacheMgr(other._cacheMgr, this), _curNormSet(0), m("m", this, other.m),
146 _varList("varList", this, other._varList), _pdfList("pdfList", this, other._pdfList), _setting(other._setting),
147 _useHorizMorph(other._useHorizMorph)
148{
149 _mref = new TVectorD(*other._mref);
152
153 // initialization
154 initialize();
155}
156
157//_____________________________________________________________________________
159{
160 if (_mref)
161 delete _mref;
162 if (_varItr)
163 delete _varItr;
164 if (_pdfItr)
165 delete _pdfItr;
166 if (_M)
167 delete _M;
168}
169
170//_____________________________________________________________________________
172{
173
174 Int_t nPdf = _pdfList.getSize();
175
176 // other quantities needed
177 if (nPdf != _mref->GetNrows()) {
178 coutE(InputArguments) << "RooMomentMorphFunc::initialize(" << GetName() << ") ERROR: nPdf != nRefPoints" << endl;
179 assert(0);
180 }
181
182 TVectorD *dm = new TVectorD(nPdf);
183 _M = new TMatrixD(nPdf, nPdf);
184
185 // transformation matrix for non-linear extrapolation, needed in evaluate()
186 TMatrixD M(nPdf, nPdf);
187 for (Int_t i = 0; i < _mref->GetNrows(); ++i) {
188 (*dm)[i] = (*_mref)[i] - (*_mref)[0];
189 M(i, 0) = 1.;
190 if (i > 0)
191 M(0, i) = 0.;
192 }
193 for (Int_t i = 1; i < _mref->GetNrows(); ++i) {
194 for (Int_t j = 1; j < _mref->GetNrows(); ++j) {
195 M(i, j) = TMath::Power((*dm)[i], (double)j);
196 }
197 }
198 (*_M) = M.Invert();
199
200 delete dm;
201}
202
203//_____________________________________________________________________________
205{
206 CacheElem *cache = (CacheElem *)_cacheMgr.getObj(0, (RooArgSet *)0);
207 if (cache) {
208 return cache;
209 }
210 Int_t nVar = _varList.getSize();
211 Int_t nPdf = _pdfList.getSize();
212
213 RooAbsReal *null = 0;
214 vector<RooAbsReal *> meanrv(nPdf * nVar, null);
215 vector<RooAbsReal *> sigmarv(nPdf * nVar, null);
216 vector<RooAbsReal *> myrms(nVar, null);
217 vector<RooAbsReal *> mypos(nVar, null);
218 vector<RooAbsReal *> slope(nPdf * nVar, null);
219 vector<RooAbsReal *> offs(nPdf * nVar, null);
220 vector<RooAbsReal *> transVar(nPdf * nVar, null);
221 vector<RooAbsReal *> transPdf(nPdf, null);
222
223 RooArgSet ownedComps;
224
225 RooArgList fracl;
226
227 // fraction parameters
228 RooArgList coefList("coefList");
229 RooArgList coefList2("coefList2");
230 for (Int_t i = 0; i < 2 * nPdf; ++i) {
231 std::string fracName = Form("frac_%d", i);
232
233 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), 1.);
234
235 fracl.add(*frac); // to be set later
236 if (i < nPdf)
237 coefList.add(*(RooRealVar *)(fracl.at(i)));
238 else
239 coefList2.add(*(RooRealVar *)(fracl.at(i)));
240 ownedComps.add(*(RooRealVar *)(fracl.at(i)));
241 }
242
243 RooRealSumFunc *theSumFunc = 0;
244 std::string sumfuncName = Form("%s_sumfunc", GetName());
245
246 if (_useHorizMorph) {
247 // mean and sigma
248 RooArgList varList(_varList);
249 for (Int_t i = 0; i < nPdf; ++i) {
250 for (Int_t j = 0; j < nVar; ++j) {
251
252 std::string meanName = Form("%s_mean_%d_%d", GetName(), i, j);
253 std::string sigmaName = Form("%s_sigma_%d_%d", GetName(), i, j);
254
255 RooAbsMoment *mom = nVar == 1 ? ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*varList.at(j))
256 : ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*varList.at(j), varList);
257
260
261 sigmarv[ij(i, j)] = mom;
262 meanrv[ij(i, j)] = mom->mean();
263
264 ownedComps.add(*sigmarv[ij(i, j)]);
265 }
266 }
267
268 // slope and offset (to be set later, depend on m)
269 for (Int_t j = 0; j < nVar; ++j) {
270 RooArgList meanList("meanList");
271 RooArgList rmsList("rmsList");
272 for (Int_t i = 0; i < nPdf; ++i) {
273 meanList.add(*meanrv[ij(i, j)]);
274 rmsList.add(*sigmarv[ij(i, j)]);
275 }
276 std::string myrmsName = Form("%s_rms_%d", GetName(), j);
277 std::string myposName = Form("%s_pos_%d", GetName(), j);
278 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList2);
279 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
280 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
281 }
282
283 // construction of unit pdfs
284 _pdfItr->Reset();
285 RooAbsPdf *pdf;
286 RooArgList transPdfList;
287
288 for (Int_t i = 0; i < nPdf; ++i) {
289 _varItr->Reset();
290 RooRealVar *var;
291
292 pdf = (RooAbsPdf *)_pdfItr->Next();
293 std::string pdfName = Form("pdf_%d", i);
294 RooCustomizer cust(*pdf, pdfName.c_str());
295
296 for (Int_t j = 0; j < nVar; ++j) {
297 // slope and offset formulas
298 std::string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
299 std::string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
300 slope[ij(i, j)] = new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[ij(i, j)], *myrms[j]));
301 offs[ij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
302 RooArgList(*meanrv[ij(i, j)], *mypos[j], *slope[ij(i, j)]));
303 ownedComps.add(RooArgSet(*slope[ij(i, j)], *offs[ij(i, j)]));
304 // linear transformations, so pdf can be renormalized
305 var = (RooRealVar *)(_varItr->Next());
306 std::string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
307 // transVar[ij(i,j)] = new
308 // RooFormulaVar(transVarName.c_str(),transVarName.c_str(),"@0*@1+@2",RooArgList(*var,*slope[ij(i,j)],*offs[ij(i,j)]));
309
310 transVar[ij(i, j)] =
311 new RooLinearVar(transVarName.c_str(), transVarName.c_str(), *var, *slope[ij(i, j)], *offs[ij(i, j)]);
312
313 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
314 // This will prevent the likelihood optimizers from erroneously declaring terms constant
315 transVar[ij(i, j)]->addServer((RooAbsArg &)m.arg());
316
317 ownedComps.add(*transVar[ij(i, j)]);
318 cust.replaceArg(*var, *transVar[ij(i, j)]);
319 }
320 transPdf[i] = (RooAbsPdf *)cust.build();
321 transPdfList.add(*transPdf[i]);
322 ownedComps.add(*transPdf[i]);
323 }
324 // sum pdf
325 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), transPdfList, coefList);
326 } else {
327 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), _pdfList, coefList);
328 }
329
330 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
331 // This will prevent the likelihood optimizers from erroneously declaring terms constant
332 theSumFunc->addServer((RooAbsArg &)m.arg());
333 theSumFunc->addOwnedComponents(ownedComps);
334
335 // change tracker for fraction parameters
336 std::string trackerName = Form("%s_frac_tracker", GetName());
337 RooChangeTracker *tracker = new RooChangeTracker(trackerName.c_str(), trackerName.c_str(), m.arg(), kTRUE);
338
339 // Store it in the cache
340 cache = new CacheElem(*theSumFunc, *tracker, fracl);
341 _cacheMgr.setObj(0, 0, cache, 0);
342
343 return cache;
344}
345
346//_____________________________________________________________________________
348{
349 return RooArgList(*_sumFunc, *_tracker);
350}
351
352//_____________________________________________________________________________
354{
355 delete _sumFunc;
356 delete _tracker;
357}
358
359//_____________________________________________________________________________
361{
362 // Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
363 _curNormSet = set ? (RooArgSet *)set : (RooArgSet *)&_varList;
364 return RooAbsReal::getVal(set);
365}
366
367//_____________________________________________________________________________
369{
370 CacheElem *cache = getCache(nset ? nset : _curNormSet);
371
372 if (cache->_tracker->hasChanged(kTRUE)) {
373 cache->calculateFractions(*this, kFALSE); // verbose turned off
374 }
375
376 return cache->_sumFunc;
377}
378
379//_____________________________________________________________________________
381{
382 CacheElem *cache = getCache(nset ? nset : _curNormSet);
383
384 if (cache->_tracker->hasChanged(kTRUE)) {
385 cache->calculateFractions(*this, kFALSE); // verbose turned off
386 }
387
388 return cache->_sumFunc;
389}
390
391//_____________________________________________________________________________
393{
395
396 if (cache->_tracker->hasChanged(kTRUE)) {
397 cache->calculateFractions(*this, kFALSE); // verbose turned off
398 }
399
400 Double_t ret = cache->_sumFunc->getVal(_pdfList.nset());
401 return ret;
402}
403
404//_____________________________________________________________________________
406{
407 return (RooRealVar *)(_frac.at(i));
408}
409
410//_____________________________________________________________________________
412{
413 return (RooRealVar *)(_frac.at(i));
414}
415
416//_____________________________________________________________________________
418{
419 Int_t nPdf = self._pdfList.getSize();
420
421 Double_t dm = self.m - (*self._mref)[0];
422
423 // fully non-linear
424 double sumposfrac = 0.;
425 for (Int_t i = 0; i < nPdf; ++i) {
426 double ffrac = 0.;
427 for (Int_t j = 0; j < nPdf; ++j) {
428 ffrac += (*self._M)(j, i) * (j == 0 ? 1. : TMath::Power(dm, (double)j));
429 }
430 if (ffrac >= 0)
431 sumposfrac += ffrac;
432 // fractions for pdf
433 ((RooRealVar *)frac(i))->setVal(ffrac);
434 // fractions for rms and mean
435 ((RooRealVar *)frac(nPdf + i))->setVal(ffrac);
436 if (verbose) {
437 cout << ffrac << endl;
438 }
439 }
440
441 // various mode settings
442 int imin = self.idxmin(self.m);
443 int imax = self.idxmax(self.m);
444 double mfrac = (self.m - (*self._mref)[imin]) / ((*self._mref)[imax] - (*self._mref)[imin]);
445 switch (self._setting) {
446 case NonLinear:
447 // default already set above
448 break;
449
450 case SineLinear:
451 mfrac =
452 TMath::Sin(TMath::PiOver2() * mfrac); // this gives a continuous differentiable transition between grid points.
453
454 // now fall through to Linear case
455
456 case Linear:
457 for (Int_t i = 0; i < 2 * nPdf; ++i) ((RooRealVar *)frac(i))->setVal(0.);
458 if (imax > imin) { // m in between mmin and mmax
459 ((RooRealVar *)frac(imin))->setVal(1. - mfrac);
460 ((RooRealVar *)frac(nPdf + imin))->setVal(1. - mfrac);
461 ((RooRealVar *)frac(imax))->setVal(mfrac);
462 ((RooRealVar *)frac(nPdf + imax))->setVal(mfrac);
463 } else if (imax == imin) { // m outside mmin and mmax
464 ((RooRealVar *)frac(imin))->setVal(1.);
465 ((RooRealVar *)frac(nPdf + imin))->setVal(1.);
466 }
467 break;
469 for (Int_t i = 0; i < nPdf; ++i) ((RooRealVar *)frac(i))->setVal(0.);
470 if (imax > imin) { // m in between mmin and mmax
471 ((RooRealVar *)frac(imin))->setVal(1. - mfrac);
472 ((RooRealVar *)frac(imax))->setVal(mfrac);
473 } else if (imax == imin) { // m outside mmin and mmax
474 ((RooRealVar *)frac(imin))->setVal(1.);
475 }
476 break;
478 for (Int_t i = 0; i < nPdf; ++i) {
479 if (((RooRealVar *)frac(i))->getVal() < 0)
480 ((RooRealVar *)frac(i))->setVal(0.);
481 ((RooRealVar *)frac(i))->setVal(((RooRealVar *)frac(i))->getVal() / sumposfrac);
482 }
483 break;
484 }
485}
486
487//_____________________________________________________________________________
488int RooMomentMorphFunc::idxmin(const double &mval) const
489{
490 int imin(0);
491 Int_t nPdf = _pdfList.getSize();
492 double mmin = -DBL_MAX;
493 for (Int_t i = 0; i < nPdf; ++i)
494 if ((*_mref)[i] > mmin && (*_mref)[i] <= mval) {
495 mmin = (*_mref)[i];
496 imin = i;
497 }
498 return imin;
499}
500
501//_____________________________________________________________________________
502int RooMomentMorphFunc::idxmax(const double &mval) const
503{
504 int imax(0);
505 Int_t nPdf = _pdfList.getSize();
506 double mmax = DBL_MAX;
507 for (Int_t i = 0; i < nPdf; ++i)
508 if ((*_mref)[i] < mmax && (*_mref)[i] >= mval) {
509 mmax = (*_mref)[i];
510 imax = i;
511 }
512 return imax;
513}
514
515//_____________________________________________________________________________
517{
518 return sumFunc(0)->plotSamplingHint(obs, xlo, xhi);
519}
520
521//_____________________________________________________________________________
522std::list<Double_t> *RooMomentMorphFunc::binBoundaries(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
523{
524 return sumFunc(0)->binBoundaries(obs, xlo, xhi);
525}
526
527//_____________________________________________________________________________
529{
530 return sumFunc(0)->isBinnedDistribution(obs);
531}
#define coutW(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#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:72
friend class RooArgSet
Definition RooAbsArg.h:606
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_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
void setLocalNoDirtyInhibit(Bool_t flag) const
Definition RooAbsArg.h:678
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
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
friend class RooRealSumFunc
Definition RooAbsReal.h:548
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:341
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:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:70
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
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 RooMomentMorphFunc &self, Bool_t verbose=kTRUE) const
virtual RooArgList containedArgs(Action)
TIterator * _pdfItr
do not persist
virtual Double_t getVal(const RooArgSet *set=0) const
RooArgSet * _curNormSet
The cache manager.
RooObjCacheManager _cacheMgr
CacheElem * getCache(const RooArgSet *nset) const
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Retrieve bin boundaries if this distribution is binned in obs.
RooAbsReal * sumFunc(const RooArgSet *nset)
Bool_t isBinnedDistribution(const RooArgSet &obs) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
int idxmin(const double &m) const
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t ij(const Int_t &i, const Int_t &j) const
int idxmax(const double &m) const
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
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:735
Double_t Sin(Double_t)
Definition TMath.h:639
auto * m
Definition textangle.C:8