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 : _cacheMgr(this, 10, true, true), _curNormSet(0), _mref(0), _M(0), _useHorizMorph(true)
34{
35}
36
37//_____________________________________________________________________________
38RooMomentMorphFunc::RooMomentMorphFunc(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
39 const RooArgList &pdfList, const TVectorD &mrefpoints, Setting setting)
40 : RooAbsReal(name, title), _cacheMgr(this, 10, true, true), m("m", "m", this, _m),
41 _varList("varList", "List of variables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
42 _useHorizMorph(true)
43{
44 // CTOR
45
46 // observables
47 for (auto *var : varList) {
48 if (!dynamic_cast<RooAbsReal *>(var)) {
49 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
50 << " is not of type RooAbsReal" << endl;
51 throw string("RooPolyMorh::ctor() ERROR variable is not of type RooAbsReal");
52 }
53 _varList.add(*var);
54 }
55
56 // reference p.d.f.s
57 for (auto const *pdf : pdfList) {
58 if (!dynamic_cast<RooAbsReal const*>(pdf)) {
59 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: func " << pdf->GetName()
60 << " is not of type RooAbsReal" << endl;
61 throw string("RooMomentMorhFunc::ctor() ERROR func is not of type RooAbsReal");
62 }
63 _pdfList.add(*pdf);
64 }
65
66 _mref = new TVectorD(mrefpoints);
67
68 // initialization
69 initialize();
70}
71
72//_____________________________________________________________________________
73RooMomentMorphFunc::RooMomentMorphFunc(const char *name, const char *title, RooAbsReal &_m, const RooArgList &varList,
74 const RooArgList &pdfList, const RooArgList &mrefList, Setting setting)
75 : RooAbsReal(name, title), _cacheMgr(this, 10, true, true), m("m", "m", this, _m),
76 _varList("varList", "List of variables", this), _pdfList("pdfList", "List of pdfs", this), _setting(setting),
77 _useHorizMorph(true)
78{
79 // CTOR
80
81 // observables
82 for (auto *var : varList) {
83 if (!dynamic_cast<RooAbsReal *>(var)) {
84 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
85 << " is not of type RooAbsReal" << endl;
86 throw string("RooMomentMorh::ctor() ERROR variable is not of type RooAbsReal");
87 }
88 _varList.add(*var);
89 }
90
91 // reference p.d.f.s
92 for (auto const *pdf : pdfList){
93 if (!dynamic_cast<RooAbsPdf const*>(pdf)) {
94 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: function " << pdf->GetName()
95 << " is not of type RooAbsReal" << endl;
96 throw string("RooMomentMorh::ctor() ERROR function is not of type RooAbsReal");
97 }
98 _pdfList.add(*pdf);
99 }
100
101 // reference points in m
102 _mref = new TVectorD(mrefList.getSize());
103 Int_t i = 0;
104 for (auto *mref : mrefList) {
105 if (!dynamic_cast<RooAbsReal *>(mref)) {
106 coutE(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") ERROR: mref " << mref->GetName()
107 << " is not of type RooAbsReal" << endl;
108 throw string("RooPolyMorh::ctor() ERROR mref is not of type RooAbsReal");
109 }
110 if (!dynamic_cast<RooConstVar *>(mref)) {
111 coutW(InputArguments) << "RooMomentMorphFunc::ctor(" << GetName() << ") WARNING mref point " << i
112 << " is not a constant, taking a snapshot of its value" << endl;
113 }
114 (*_mref)[i] = static_cast<RooAbsReal *>(mref)->getVal();
115 ++i;
116 }
117
118 // initialization
119 initialize();
120}
121
122//_____________________________________________________________________________
124 : RooAbsReal(other, name), _cacheMgr(other._cacheMgr, this), _curNormSet(0), m("m", this, other.m),
125 _varList("varList", this, other._varList), _pdfList("pdfList", this, other._pdfList), _setting(other._setting),
126 _useHorizMorph(other._useHorizMorph)
127{
128 _mref = new TVectorD(*other._mref);
129
130 // initialization
131 initialize();
132}
133
134//_____________________________________________________________________________
136{
137 if (_mref)
138 delete _mref;
139 if (_M)
140 delete _M;
141}
142
143//_____________________________________________________________________________
145{
146
147 Int_t nPdf = _pdfList.getSize();
148
149 // other quantities needed
150 if (nPdf != _mref->GetNrows()) {
151 coutE(InputArguments) << "RooMomentMorphFunc::initialize(" << GetName() << ") ERROR: nPdf != nRefPoints" << endl;
152 assert(0);
153 }
154
155 TVectorD *dm = new TVectorD(nPdf);
156 _M = new TMatrixD(nPdf, nPdf);
157
158 // transformation matrix for non-linear extrapolation, needed in evaluate()
159 TMatrixD M(nPdf, nPdf);
160 for (Int_t i = 0; i < _mref->GetNrows(); ++i) {
161 (*dm)[i] = (*_mref)[i] - (*_mref)[0];
162 M(i, 0) = 1.;
163 if (i > 0)
164 M(0, i) = 0.;
165 }
166 for (Int_t i = 1; i < _mref->GetNrows(); ++i) {
167 for (Int_t j = 1; j < _mref->GetNrows(); ++j) {
168 M(i, j) = TMath::Power((*dm)[i], (double)j);
169 }
170 }
171 (*_M) = M.Invert();
172
173 delete dm;
174}
175
176//_____________________________________________________________________________
178{
179 CacheElem *cache = (CacheElem *)_cacheMgr.getObj(0, (RooArgSet *)0);
180 if (cache) {
181 return cache;
182 }
183 Int_t nVar = _varList.getSize();
184 Int_t nPdf = _pdfList.getSize();
185
186 RooAbsReal *null = 0;
187 vector<RooAbsReal *> meanrv(nPdf * nVar, null);
188 vector<RooAbsReal *> sigmarv(nPdf * nVar, null);
189 vector<RooAbsReal *> myrms(nVar, null);
190 vector<RooAbsReal *> mypos(nVar, null);
191 vector<RooAbsReal *> slope(nPdf * nVar, null);
192 vector<RooAbsReal *> offs(nPdf * nVar, null);
193 vector<RooAbsReal *> transVar(nPdf * nVar, null);
194 vector<RooAbsReal *> transPdf(nPdf, null);
195
196 RooArgSet ownedComps;
197
198 RooArgList fracl;
199
200 // fraction parameters
201 RooArgList coefList("coefList");
202 RooArgList coefList2("coefList2");
203 for (Int_t i = 0; i < 2 * nPdf; ++i) {
204 std::string fracName = Form("frac_%d", i);
205
206 RooRealVar *frac = new RooRealVar(fracName.c_str(), fracName.c_str(), 1.);
207
208 fracl.add(*frac); // to be set later
209 if (i < nPdf)
210 coefList.add(*(RooRealVar *)(fracl.at(i)));
211 else
212 coefList2.add(*(RooRealVar *)(fracl.at(i)));
213 ownedComps.add(*(RooRealVar *)(fracl.at(i)));
214 }
215
216 RooRealSumFunc *theSumFunc = 0;
217 std::string sumfuncName = Form("%s_sumfunc", GetName());
218
219 if (_useHorizMorph) {
220 // mean and sigma
221 RooArgList varList(_varList);
222 for (Int_t i = 0; i < nPdf; ++i) {
223 for (Int_t j = 0; j < nVar; ++j) {
224
225 std::string meanName = Form("%s_mean_%d_%d", GetName(), i, j);
226 std::string sigmaName = Form("%s_sigma_%d_%d", GetName(), i, j);
227
228 RooAbsMoment *mom = nVar == 1 ? ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*varList.at(j))
229 : ((RooAbsPdf *)_pdfList.at(i))->sigma((RooRealVar &)*varList.at(j), varList);
230
231 mom->setLocalNoDirtyInhibit(true);
232 mom->mean()->setLocalNoDirtyInhibit(true);
233
234 sigmarv[ij(i, j)] = mom;
235 meanrv[ij(i, j)] = mom->mean();
236
237 ownedComps.add(*sigmarv[ij(i, j)]);
238 }
239 }
240
241 // slope and offset (to be set later, depend on m)
242 for (Int_t j = 0; j < nVar; ++j) {
243 RooArgList meanList("meanList");
244 RooArgList rmsList("rmsList");
245 for (Int_t i = 0; i < nPdf; ++i) {
246 meanList.add(*meanrv[ij(i, j)]);
247 rmsList.add(*sigmarv[ij(i, j)]);
248 }
249 std::string myrmsName = Form("%s_rms_%d", GetName(), j);
250 std::string myposName = Form("%s_pos_%d", GetName(), j);
251 myrms[j] = new RooAddition(myrmsName.c_str(), myrmsName.c_str(), rmsList, coefList2);
252 mypos[j] = new RooAddition(myposName.c_str(), myposName.c_str(), meanList, coefList2);
253 ownedComps.add(RooArgSet(*myrms[j], *mypos[j]));
254 }
255
256 // construction of unit pdfs
257 RooArgList transPdfList;
258
259 for (Int_t i = 0; i < nPdf; ++i) {
260 auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
261 std::string pdfName = Form("pdf_%d", i);
262 RooCustomizer cust(pdf, pdfName.c_str());
263
264 for (Int_t j = 0; j < nVar; ++j) {
265 // slope and offset formulas
266 std::string slopeName = Form("%s_slope_%d_%d", GetName(), i, j);
267 std::string offsetName = Form("%s_offset_%d_%d", GetName(), i, j);
268 slope[ij(i, j)] = new RooFormulaVar(slopeName.c_str(), "@0/@1", RooArgList(*sigmarv[ij(i, j)], *myrms[j]));
269 offs[ij(i, j)] = new RooFormulaVar(offsetName.c_str(), "@0-(@1*@2)",
270 RooArgList(*meanrv[ij(i, j)], *mypos[j], *slope[ij(i, j)]));
271 ownedComps.add(RooArgSet(*slope[ij(i, j)], *offs[ij(i, j)]));
272 // linear transformations, so pdf can be renormalized
273 auto& var = static_cast<RooRealVar&>(*_varList[j]);
274 std::string transVarName = Form("%s_transVar_%d_%d", GetName(), i, j);
275 // transVar[ij(i,j)] = new
276 // RooFormulaVar(transVarName.c_str(),transVarName.c_str(),"@0*@1+@2",RooArgList(*var,*slope[ij(i,j)],*offs[ij(i,j)]));
277
278 transVar[ij(i, j)] =
279 new RooLinearVar(transVarName.c_str(), transVarName.c_str(), var, *slope[ij(i, j)], *offs[ij(i, j)]);
280
281 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
282 // This will prevent the likelihood optimizers from erroneously declaring terms constant
283 transVar[ij(i, j)]->addServer((RooAbsArg &)m.arg());
284
285 ownedComps.add(*transVar[ij(i, j)]);
286 cust.replaceArg(var, *transVar[ij(i, j)]);
287 }
288 transPdf[i] = (RooAbsPdf *)cust.build();
289 transPdfList.add(*transPdf[i]);
290 ownedComps.add(*transPdf[i]);
291 }
292 // sum pdf
293 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), transPdfList, coefList);
294 } else {
295 theSumFunc = new RooRealSumFunc(sumfuncName.c_str(), sumfuncName.c_str(), _pdfList, coefList);
296 }
297
298 // *** WVE this is important *** this declares that frac effectively depends on the morphing parameters
299 // This will prevent the likelihood optimizers from erroneously declaring terms constant
300 theSumFunc->addServer((RooAbsArg &)m.arg());
301 theSumFunc->addOwnedComponents(ownedComps);
302
303 // change tracker for fraction parameters
304 std::string trackerName = Form("%s_frac_tracker", GetName());
305 RooChangeTracker *tracker = new RooChangeTracker(trackerName.c_str(), trackerName.c_str(), m.arg(), true);
306
307 // Store it in the cache
308 cache = new CacheElem(*theSumFunc, *tracker, fracl);
309 _cacheMgr.setObj(0, 0, cache, 0);
310
311 return cache;
312}
313
314//_____________________________________________________________________________
316{
317 return RooArgList(*_sumFunc, *_tracker);
318}
319
320//_____________________________________________________________________________
322{
323 delete _sumFunc;
324 delete _tracker;
325}
326
327//_____________________________________________________________________________
328double RooMomentMorphFunc::getVal(const RooArgSet *set) const
329{
330 // Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
331 _curNormSet = set ? (RooArgSet *)set : (RooArgSet *)&_varList;
332 return RooAbsReal::getVal(set);
333}
334
335//_____________________________________________________________________________
337{
338 CacheElem *cache = getCache(nset ? nset : _curNormSet);
339
340 if (cache->_tracker->hasChanged(true)) {
341 cache->calculateFractions(*this, false); // verbose turned off
342 }
343
344 return cache->_sumFunc;
345}
346
347//_____________________________________________________________________________
349{
350 CacheElem *cache = getCache(nset ? nset : _curNormSet);
351
352 if (cache->_tracker->hasChanged(true)) {
353 cache->calculateFractions(*this, false); // verbose turned off
354 }
355
356 return cache->_sumFunc;
357}
358
359//_____________________________________________________________________________
361{
363
364 if (cache->_tracker->hasChanged(true)) {
365 cache->calculateFractions(*this, false); // verbose turned off
366 }
367
368 double ret = cache->_sumFunc->getVal(_pdfList.nset());
369 return ret;
370}
371
372//_____________________________________________________________________________
374{
375 return (RooRealVar *)(_frac.at(i));
376}
377
378//_____________________________________________________________________________
380{
381 return (RooRealVar *)(_frac.at(i));
382}
383
384//_____________________________________________________________________________
386{
387 Int_t nPdf = self._pdfList.getSize();
388
389 double dm = self.m - (*self._mref)[0];
390
391 // fully non-linear
392 double sumposfrac = 0.;
393 for (Int_t i = 0; i < nPdf; ++i) {
394 double ffrac = 0.;
395 for (Int_t j = 0; j < nPdf; ++j) {
396 ffrac += (*self._M)(j, i) * (j == 0 ? 1. : TMath::Power(dm, (double)j));
397 }
398 if (ffrac >= 0)
399 sumposfrac += ffrac;
400 // fractions for pdf
401 ((RooRealVar *)frac(i))->setVal(ffrac);
402 // fractions for rms and mean
403 ((RooRealVar *)frac(nPdf + i))->setVal(ffrac);
404 if (verbose) {
405 cout << ffrac << endl;
406 }
407 }
408
409 // various mode settings
410 int imin = self.idxmin(self.m);
411 int imax = self.idxmax(self.m);
412 double mfrac = (self.m - (*self._mref)[imin]) / ((*self._mref)[imax] - (*self._mref)[imin]);
413 switch (self._setting) {
414 case NonLinear:
415 // default already set above
416 break;
417
418 case SineLinear:
419 mfrac =
420 TMath::Sin(TMath::PiOver2() * mfrac); // this gives a continuous differentiable transition between grid points.
421
422 // now fall through to Linear case
423
424 case Linear:
425 for (Int_t i = 0; i < 2 * nPdf; ++i) ((RooRealVar *)frac(i))->setVal(0.);
426 if (imax > imin) { // m in between mmin and mmax
427 ((RooRealVar *)frac(imin))->setVal(1. - mfrac);
428 ((RooRealVar *)frac(nPdf + imin))->setVal(1. - mfrac);
429 ((RooRealVar *)frac(imax))->setVal(mfrac);
430 ((RooRealVar *)frac(nPdf + imax))->setVal(mfrac);
431 } else if (imax == imin) { // m outside mmin and mmax
432 ((RooRealVar *)frac(imin))->setVal(1.);
433 ((RooRealVar *)frac(nPdf + imin))->setVal(1.);
434 }
435 break;
437 for (Int_t i = 0; i < nPdf; ++i) ((RooRealVar *)frac(i))->setVal(0.);
438 if (imax > imin) { // m in between mmin and mmax
439 ((RooRealVar *)frac(imin))->setVal(1. - mfrac);
440 ((RooRealVar *)frac(imax))->setVal(mfrac);
441 } else if (imax == imin) { // m outside mmin and mmax
442 ((RooRealVar *)frac(imin))->setVal(1.);
443 }
444 break;
446 for (Int_t i = 0; i < nPdf; ++i) {
447 if (((RooRealVar *)frac(i))->getVal() < 0)
448 ((RooRealVar *)frac(i))->setVal(0.);
449 ((RooRealVar *)frac(i))->setVal(((RooRealVar *)frac(i))->getVal() / sumposfrac);
450 }
451 break;
452 }
453}
454
455//_____________________________________________________________________________
456int RooMomentMorphFunc::idxmin(const double &mval) const
457{
458 int imin(0);
459 Int_t nPdf = _pdfList.getSize();
460 double mmin = -DBL_MAX;
461 for (Int_t i = 0; i < nPdf; ++i)
462 if ((*_mref)[i] > mmin && (*_mref)[i] <= mval) {
463 mmin = (*_mref)[i];
464 imin = i;
465 }
466 return imin;
467}
468
469//_____________________________________________________________________________
470int RooMomentMorphFunc::idxmax(const double &mval) const
471{
472 int imax(0);
473 Int_t nPdf = _pdfList.getSize();
474 double mmax = DBL_MAX;
475 for (Int_t i = 0; i < nPdf; ++i)
476 if ((*_mref)[i] < mmax && (*_mref)[i] >= mval) {
477 mmax = (*_mref)[i];
478 imax = i;
479 }
480 return imax;
481}
482
483//_____________________________________________________________________________
484std::list<double> *RooMomentMorphFunc::plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
485{
486 return sumFunc(0)->plotSamplingHint(obs, xlo, xhi);
487}
488
489//_____________________________________________________________________________
490std::list<double> *RooMomentMorphFunc::binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
491{
492 return sumFunc(0)->binBoundaries(obs, xlo, xhi);
493}
494
495//_____________________________________________________________________________
497{
498 return sumFunc(0)->isBinnedDistribution(obs);
499}
#define coutW(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
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:2467
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:74
void setLocalNoDirtyInhibit(bool flag) const
Definition RooAbsArg.h:697
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
void addServer(RooAbsArg &server, bool valueProp=true, bool shapeProp=false, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsMoment represents the first, second, or third order derivative of any RooAbsReal as calculated ...
RooAbsReal * mean()
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
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:62
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:91
friend class RooRealSumFunc
Definition RooAbsReal.h:492
virtual bool isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:358
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...
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:55
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.
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
bool hasChanged(bool clearState)
Returns true if state has changed since last call with clearState=true.
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...
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 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 RooMomentMorphFunc &self, bool verbose=true) const
RooArgList containedArgs(Action) override
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...
virtual double getVal(const RooArgSet *set=nullptr) const
Int_t ij(const Int_t &i, const Int_t &j) const
int idxmax(const double &m) const
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
const T & arg() const
Return reference to object held in proxy.
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
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
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)
Returns x raised to the power y.
Definition TMath.h:719
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
TMarker m
Definition textangle.C:8