Logo ROOT  
Reference Guide
RooRealSumPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/** \class RooRealSumPdf
19 \ingroup Roofitcore
20
21
22The class RooRealSumPdf implements a PDF constructed from a sum of functions:
23\f[
24 \mathrm{PDF}(x) = \frac{ \sum_{i=1}^{n-1} \mathrm{coef}_i * \mathrm{func}_i(x) + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \mathrm{func}_n(x) }
25 {\sum_{i=1}^{n-1} \mathrm{coef}_i * \int \mathrm{func}_i(x)dx + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \int \mathrm{func}_n(x) dx }
26\f]
27
28where \f$\mathrm{coef}_i\f$ and \f$\mathrm{func}_i\f$ are RooAbsReal objects, and \f$ x \f$ is the collection of dependents.
29In the present version \f$\mathrm{coef}_i\f$ may not depend on \f$ x \f$, but this limitation could be removed should the need arise.
30
31If the number of coefficients is one less than the number of functions, the PDF is assumed to be normalised. Due to this additional constraint,
32\f$\mathrm{coef}_n\f$ is computed from the other coefficients.
33
34### Extending the PDF
35If an \f$ n^\mathrm{th} \f$ coefficient is provided, the PDF **can** be used as an extended PDF, *i.e.* the total number of events will be measured in addition
36to the fractions of the various functions. **This requires setting the last argument of the constructor to `true`.**
37\note For the RooAddPdf, the extension happens automatically.
38
39### Difference to RooAddPdf / RooRealSumFunc
40- RooAddPdf is a PDF of PDFs, *i.e.* its components need to be normalised and non-negative.
41- RooRealSumPdf is a PDF of functions, *i.e.*, its components can be negative, but their sum cannot be. The normalisation
42 is computed automatically, unless the PDF is extended (see above).
43- RooRealSumFunc is a sum of functions. It is neither normalised, nor need it be positive.
44
45*/
46
47#include "RooFit.h"
48#include "Riostream.h"
49
50#include "TError.h"
51#include "TIterator.h"
52#include "TList.h"
53#include "TClass.h"
54#include "RooRealSumPdf.h"
55#include "RooRealProxy.h"
56#include "RooPlot.h"
57#include "RooRealVar.h"
58#include "RooAddGenContext.h"
59#include "RooRealConstant.h"
60#include "RooRealIntegral.h"
61#include "RooMsgService.h"
62#include "RooNameReg.h"
63
64#include <algorithm>
65#include <memory>
66
67using namespace std;
68
70;
71
73
74////////////////////////////////////////////////////////////////////////////////
75/// Default constructor
76/// coverity[UNINIT_CTOR]
77
79{
82}
83
84
85
86////////////////////////////////////////////////////////////////////////////////
87/// Constructor with name and title
88
89RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
90 RooAbsPdf(name,title),
91 _normIntMgr(this,10),
92 _funcList("!funcList","List of functions",this),
93 _coefList("!coefList","List of coefficients",this),
94 _extended(kFALSE),
95 _doFloor(kFALSE)
96{
97
98}
99
100
101
102////////////////////////////////////////////////////////////////////////////////
103/// Construct p.d.f consisting of \f$ \mathrm{coef}_1 * \mathrm{func}_1 + (1-\mathrm{coef}_1) * \mathrm{func}_2 \f$.
104/// The input coefficients and functions are allowed to be negative
105/// but the resulting sum is not, which is enforced at runtime.
106
107RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
108 RooAbsReal& func1, RooAbsReal& func2, RooAbsReal& coef1) :
109 RooAbsPdf(name,title),
110 _normIntMgr(this,10),
111 _funcList("!funcList","List of functions",this),
112 _coefList("!coefList","List of coefficients",this),
113 _extended(kFALSE),
114 _doFloor(kFALSE)
115{
116 // Special constructor with two functions and one coefficient
117
118 _funcList.add(func1) ;
119 _funcList.add(func2) ;
120 _coefList.add(coef1) ;
121
122}
123
124
125////////////////////////////////////////////////////////////////////////////////
126/// Constructor for a PDF from a list of functions and coefficients.
127/// It implements
128/// \f[
129/// \sum_i \mathrm{coef}_i \cdot \mathrm{func}_i,
130/// \f]
131/// if \f$ N_\mathrm{coef} = N_\mathrm{func} \f$. With `extended=true`, the coefficients can take any values. With `extended=false`,
132/// there is the danger of getting a degenerate minimisation problem because a PDF has to be normalised, which needs one degree
133/// of freedom less.
134///
135/// A plain (normalised) PDF can therefore be implemented with one less coefficient. RooFit then computes
136/// \f[
137/// \sum_i^{N-1} \mathrm{coef}_i \cdot \mathrm{func}_i + (1 - \sum_i \mathrm{coef}_i ) \cdot \mathrm{func}_N,
138/// \f]
139/// if \f$ N_\mathrm{coef} = N_\mathrm{func} - 1 \f$.
140///
141/// All coefficients and functions are allowed to be negative
142/// but the sum (*i.e.* the PDF) is not, which is enforced at runtime.
143///
144/// \param name Name of the PDF
145/// \param title Title (for plotting)
146/// \param inFuncList List of functions to sum
147/// \param inCoefList List of coefficients
148/// \param extended Interpret as extended PDF (requires equal number of functions and coefficients)
149
150RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
151 const RooArgList& inFuncList, const RooArgList& inCoefList, Bool_t extended) :
152 RooAbsPdf(name,title),
153 _normIntMgr(this,10),
154 _funcList("!funcList","List of functions",this),
155 _coefList("!coefList","List of coefficients",this),
156 _extended(extended),
157 _doFloor(kFALSE)
158{
159 if (!(inFuncList.getSize()==inCoefList.getSize()+1 || inFuncList.getSize()==inCoefList.getSize())) {
160 coutE(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName()
161 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << endl ;
162 assert(0) ;
163 }
164
165 // Constructor with N functions and N or N-1 coefs
166 for (unsigned int i = 0; i < inCoefList.size(); ++i) {
167 const auto& func = inFuncList[i];
168 const auto& coef = inCoefList[i];
169
170 if (!dynamic_cast<const RooAbsReal*>(&coef)) {
171 coutW(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") coefficient " << coef.GetName() << " is not of type RooAbsReal, ignored" << endl ;
172 continue ;
173 }
174 if (!dynamic_cast<const RooAbsReal*>(&func)) {
175 coutW(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") func " << func.GetName() << " is not of type RooAbsReal, ignored" << endl ;
176 continue ;
177 }
178 _funcList.add(func) ;
179 _coefList.add(coef) ;
180 }
181
182 if (inFuncList.size() == inCoefList.size() + 1) {
183 const auto& func = inFuncList[inFuncList.size()-1];
184 if (!dynamic_cast<const RooAbsReal*>(&func)) {
185 coutE(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") last func " << func.GetName() << " is not of type RooAbsReal, fatal error" << endl ;
186 assert(0) ;
187 }
188 _funcList.add(func);
189 }
190
191}
192
193
194
195
196////////////////////////////////////////////////////////////////////////////////
197/// Copy constructor
198
200 RooAbsPdf(other,name),
201 _normIntMgr(other._normIntMgr,this),
202 _funcList("!funcList",this,other._funcList),
203 _coefList("!coefList",this,other._coefList),
204 _extended(other._extended),
205 _doFloor(other._doFloor)
206{
207
208}
209
210
211
212////////////////////////////////////////////////////////////////////////////////
213/// Destructor
214
216{
217
218}
219
220
221
222
223
224////////////////////////////////////////////////////////////////////////////////
225
227{
229}
230
231
232
233
234////////////////////////////////////////////////////////////////////////////////
235/// Calculate the current value
236
238{
239 Double_t value(0) ;
240
241 // Do running sum of coef/func pairs, calculate lastCoef.
242
243 // N funcs, N-1 coefficients
244 Double_t lastCoef(1) ;
245 auto funcIt = _funcList.begin();
246 for (const auto coefArg : _coefList) {
247 assert(funcIt != _funcList.end());
248 auto func = static_cast<const RooAbsReal*>(*funcIt++);
249 auto coef = static_cast<const RooAbsReal*>(coefArg);
250
251 Double_t coefVal = coef->getVal() ;
252 if (coefVal) {
253 cxcoutD(Eval) << "RooRealSumPdf::eval(" << GetName() << ") coefVal = " << coefVal << " funcVal = " << func->IsA()->GetName() << "::" << func->GetName() << " = " << func->getVal() << endl ;
254 if (func->isSelectedComp()) {
255 value += func->getVal()*coefVal ;
256 }
257 lastCoef -= coef->getVal() ;
258 }
259 }
260
261 if (!haveLastCoef()) {
262 assert(funcIt != _funcList.end());
263 // Add last func with correct coefficient
264 auto func = static_cast<const RooAbsReal*>(*funcIt);
265 if (func->isSelectedComp()) {
266 value += func->getVal()*lastCoef ;
267 }
268
269 cxcoutD(Eval) << "RooRealSumPdf::eval(" << GetName() << ") lastCoef = " << lastCoef << " funcVal = " << func->getVal() << endl ;
270
271 // Warn about coefficient degeneration
272 if (lastCoef<0 || lastCoef>1) {
273 coutW(Eval) << "RooRealSumPdf::evaluate(" << GetName()
274 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
275 << 1-lastCoef << ". This means that the PDF is not properly normalised. If the PDF was meant to be extended, provide as many coefficients as functions." << endl ;
276 }
277 }
278
279 // Introduce floor if so requested
280 if (value<0 && (_doFloor || _doFloorGlobal)) {
281 value = 0 ;
282 }
283
284 return value ;
285}
286
287
288
289
290////////////////////////////////////////////////////////////////////////////////
291/// Check if FUNC is valid for given normalization set.
292/// Coefficient and FUNC must be non-overlapping, but func-coefficient
293/// pairs may overlap each other.
294///
295/// In the present implementation, coefficients may not be observables or derive
296/// from observables.
297
299{
300 Bool_t ret(kFALSE) ;
301
302 for (unsigned int i=0; i < _coefList.size(); ++i) {
303 const auto& coef = _coefList[i];
304 const auto& func = _funcList[i];
305
306 if (func.observableOverlaps(nset, coef)) {
307 coutE(InputArguments) << "RooRealSumPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef.GetName()
308 << " and FUNC " << func.GetName() << " have one or more observables in common" << endl ;
309 ret = kTRUE ;
310 }
311 if (coef.dependsOn(*nset)) {
312 coutE(InputArguments) << "RooRealPdf::checkObservables(" << GetName() << "): ERROR coefficient " << coef.GetName()
313 << " depends on one or more of the following observables" ; nset->Print("1") ;
314 ret = kTRUE ;
315 }
316 }
317
318 return ret ;
319}
320
321
322
323
324////////////////////////////////////////////////////////////////////////////////
325/// Advertise that all integrals can be handled internally.
326
328 const RooArgSet* normSet2, const char* rangeName) const
329{
330 // Handle trivial no-integration scenario
331 if (allVars.getSize()==0) return 0 ;
332 if (_forceNumInt) return 0 ;
333
334 // Select subset of allVars that are actual dependents
335 analVars.add(allVars) ;
336 RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
337
338
339 // Check if this configuration was created before
340 Int_t sterileIdx(-1) ;
341 CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx,RooNameReg::ptr(rangeName)) ;
342 if (cache) {
343 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
344 return _normIntMgr.lastIndex()+1 ;
345 }
346
347 // Create new cache element
348 cache = new CacheElem ;
349
350 // Make list of function projection and normalization integrals
351 for (const auto elm : _funcList) {
352 const auto func = static_cast<RooAbsReal*>(elm);
353
354 RooAbsReal* funcInt = func->createIntegral(analVars,rangeName) ;
355 if(funcInt->InheritsFrom(RooRealIntegral::Class())) ((RooRealIntegral*)funcInt)->setAllowComponentSelection(true);
356 cache->_funcIntList.addOwned(*funcInt) ;
357 if (normSet && normSet->getSize()>0) {
358 RooAbsReal* funcNorm = func->createIntegral(*normSet) ;
359 cache->_funcNormList.addOwned(*funcNorm) ;
360 }
361 }
362
363 // Store cache element
364 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
365
366 if (normSet) {
367 delete normSet ;
368 }
369
370 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
371 return code+1 ;
372}
373
374
375
376
377////////////////////////////////////////////////////////////////////////////////
378/// Implement analytical integrations by deferring integration of component
379/// functions to integrators of components.
380
381Double_t RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
382{
383 // Handle trivial passthrough scenario
384 if (code==0) return getVal(normSet2) ;
385
386
387 // WVE needs adaptation for rangeName feature
388 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
389 if (cache==0) { // revive the (sterilized) cache
390 //cout << "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << ": reviving cache "<< endl;
391 std::unique_ptr<RooArgSet> vars( getParameters(RooArgSet()) );
392 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars) );
393 std::unique_ptr<RooArgSet> nset( _normIntMgr.nameSet1ByIndex(code-1)->select(*vars) );
395 Int_t code2 = getAnalyticalIntegralWN(*iset,dummy,nset.get(),rangeName);
396 R__ASSERT(code==code2); // must have revived the right (sterilized) slot...
397 cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
398 R__ASSERT(cache!=0);
399 }
400
401 Double_t value(0) ;
402
403 // N funcs, N-1 coefficients
404 Double_t lastCoef(1) ;
405 auto funcIt = _funcList.begin();
406 auto funcIntIt = cache->_funcIntList.begin();
407 for (const auto coefArg : _coefList) {
408 assert(funcIt != _funcList.end());
409 const auto coef = static_cast<const RooAbsReal*>(coefArg);
410 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
411 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
412
413 Double_t coefVal = coef->getVal(normSet2) ;
414 if (coefVal) {
415 assert(func);
416 if (normSet2 ==0 || func->isSelectedComp()) {
417 assert(funcInt);
418 value += funcInt->getVal()*coefVal ;
419 }
420 lastCoef -= coef->getVal(normSet2) ;
421 }
422 }
423
424 if (!haveLastCoef()) {
425 // Add last func with correct coefficient
426 const auto func = static_cast<const RooAbsReal*>(*funcIt);
427 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
428 assert(func);
429
430 if (normSet2 ==0 || func->isSelectedComp()) {
431 assert(funcInt);
432 value += funcInt->getVal()*lastCoef ;
433 }
434
435 // Warn about coefficient degeneration
436 if (lastCoef<0 || lastCoef>1) {
437 coutW(Eval) << "RooRealSumPdf::evaluate(" << GetName()
438 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
439 << 1-lastCoef << endl ;
440 }
441 }
442
443 Double_t normVal(1) ;
444 if (normSet2 && normSet2->getSize()>0) {
445 normVal = 0 ;
446
447 // N funcs, N-1 coefficients
448 auto funcNormIter = cache->_funcNormList.begin();
449 for (const auto coefAsArg : _coefList) {
450 auto coef = static_cast<RooAbsReal*>(coefAsArg);
451 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
452
453 Double_t coefVal = coef->getVal(normSet2);
454 if (coefVal) {
455 assert(funcNorm);
456 normVal += funcNorm->getVal()*coefVal ;
457 }
458 }
459
460 // Add last func with correct coefficient
461 if (!haveLastCoef()) {
462 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
463 assert(funcNorm);
464
465 normVal += funcNorm->getVal()*lastCoef;
466 }
467 }
468
469 return value / normVal;
470}
471
472
473////////////////////////////////////////////////////////////////////////////////
474
476{
477 Double_t n = getNorm(nset) ;
478 if (n<0) {
479 logEvalError("Expected number of events is negative") ;
480 }
481 return n ;
482}
483
484
485////////////////////////////////////////////////////////////////////////////////
486
487std::list<Double_t>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
488{
489 list<Double_t>* sumBinB = 0 ;
490 Bool_t needClean(kFALSE) ;
491
492 // Loop over components pdf
493 for (const auto elm : _funcList) {
494 auto func = static_cast<RooAbsReal*>(elm);
495
496 list<Double_t>* funcBinB = func->binBoundaries(obs,xlo,xhi) ;
497
498 // Process hint
499 if (funcBinB) {
500 if (!sumBinB) {
501 // If this is the first hint, then just save it
502 sumBinB = funcBinB ;
503 } else {
504
505 list<Double_t>* newSumBinB = new list<Double_t>(sumBinB->size()+funcBinB->size()) ;
506
507 // Merge hints into temporary array
508 merge(funcBinB->begin(),funcBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
509
510 // Copy merged array without duplicates to new sumBinBArrau
511 delete sumBinB ;
512 delete funcBinB ;
513 sumBinB = newSumBinB ;
514 needClean = kTRUE ;
515 }
516 }
517 }
518
519 // Remove consecutive duplicates
520 if (needClean) {
521 list<Double_t>::iterator new_end = unique(sumBinB->begin(),sumBinB->end()) ;
522 sumBinB->erase(new_end,sumBinB->end()) ;
523 }
524
525 return sumBinB ;
526}
527
528
529
530//_____________________________________________________________________________B
532{
533 // If all components that depend on obs are binned that so is the product
534
535 for (const auto elm : _funcList) {
536 auto func = static_cast<RooAbsReal*>(elm);
537
538 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
539 return kFALSE ;
540 }
541 }
542
543 return kTRUE ;
544}
545
546
547
548
549
550////////////////////////////////////////////////////////////////////////////////
551
552std::list<Double_t>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
553{
554 list<Double_t>* sumHint = 0 ;
555 Bool_t needClean(kFALSE) ;
556
557 // Loop over components pdf
558 for (const auto elm : _funcList) {
559 auto func = static_cast<RooAbsReal*>(elm);
560
561 list<Double_t>* funcHint = func->plotSamplingHint(obs,xlo,xhi) ;
562
563 // Process hint
564 if (funcHint) {
565 if (!sumHint) {
566
567 // If this is the first hint, then just save it
568 sumHint = funcHint ;
569
570 } else {
571
572 list<Double_t>* newSumHint = new list<Double_t>(sumHint->size()+funcHint->size()) ;
573
574 // Merge hints into temporary array
575 merge(funcHint->begin(),funcHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
576
577 // Copy merged array without duplicates to new sumHintArrau
578 delete sumHint ;
579 sumHint = newSumHint ;
580 needClean = kTRUE ;
581 }
582 }
583 }
584
585 // Remove consecutive duplicates
586 if (needClean) {
587 list<Double_t>::iterator new_end = unique(sumHint->begin(),sumHint->end()) ;
588 sumHint->erase(new_end,sumHint->end()) ;
589 }
590
591 return sumHint ;
592}
593
594
595
596
597////////////////////////////////////////////////////////////////////////////////
598/// Label OK'ed components of a RooRealSumPdf with cache-and-track
599
601{
602 for (const auto sarg : _funcList) {
603 if (sarg->canNodeBeCached()==Always) {
604 trackNodes.add(*sarg) ;
605 //cout << "tracking node RealSumPdf component " << sarg->IsA()->GetName() << "::" << sarg->GetName() << endl ;
606 }
607 }
608}
609
610
611
612////////////////////////////////////////////////////////////////////////////////
613/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
614/// product operator construction
615
616void RooRealSumPdf::printMetaArgs(ostream& os) const
617{
618
620
621 if (_coefList.getSize()!=0) {
622 auto funcIter = _funcList.begin();
623
624 for (const auto coef : _coefList) {
625 if (!first) {
626 os << " + " ;
627 } else {
628 first = kFALSE ;
629 }
630 const auto func = *(funcIter++);
631 os << coef->GetName() << " * " << func->GetName();
632 }
633
634 if (funcIter != _funcList.end()) {
635 os << " + [%] * " << (*funcIter)->GetName() ;
636 }
637 } else {
638
639 for (const auto func : _funcList) {
640 if (!first) {
641 os << " + " ;
642 } else {
643 first = kFALSE ;
644 }
645 os << func->GetName() ;
646 }
647 }
648
649 os << " " ;
650}
void Class()
Definition: Class.C:29
static RooMathCoreReg dummy
#define cxcoutD(a)
Definition: RooMsgService.h:81
#define coutW(a)
Definition: RooMsgService.h:32
#define coutE(a)
Definition: RooMsgService.h:33
const Bool_t kFALSE
Definition: RtypesCore.h:90
bool Bool_t
Definition: RtypesCore.h:61
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define ClassImp(name)
Definition: Rtypes.h:361
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:276
friend class RooArgSet
Definition: RooAbsArg.h:572
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:544
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
Int_t getSize() const
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
const_iterator end() const
Storage_t::size_type size() const
const_iterator begin() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
Double_t getNorm(const RooArgSet &nset) const
Definition: RooAbsPdf.h:207
friend class CacheElem
Definition: RooAbsPdf.h:333
@ CanBeExtended
Definition: RooAbsPdf.h:229
@ CanNotBeExtended
Definition: RooAbsPdf.h:229
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:60
Bool_t _forceNumInt
Definition: RooAbsReal.h:453
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:314
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Definition: RooAbsReal.h:313
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:560
void logEvalError(const char *message, const char *serverValueString=0) const
Log evaluation error message.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
const RooNameSet * nameSet2ByIndex(Int_t index) const
Retrieve RooNameSet associated with slot at given index.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
Int_t lastIndex() const
const RooNameSet * nameSet1ByIndex(Int_t index) const
Retrieve RooNameSet associated with slot at given index.
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:93
RooArgSet * select(const RooArgSet &list) const
Construct a RooArgSet of objects in input 'list' whose names match to those in the internal name list...
Definition: RooNameSet.cxx:185
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
Definition: RooRealSumPdf.h:24
virtual ExtendMode extendMode() const
RooObjCacheManager _normIntMgr
Definition: RooRealSumPdf.h:82
Double_t evaluate() const
Calculate the current value.
RooListProxy _coefList
Definition: RooRealSumPdf.h:86
RooListProxy _funcList
Definition: RooRealSumPdf.h:85
RooRealSumPdf()
Default constructor coverity[UNINIT_CTOR].
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Bool_t isBinnedDistribution(const RooArgSet &obs) const
static Bool_t _doFloorGlobal
Definition: RooRealSumPdf.h:90
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
virtual ~RooRealSumPdf()
Destructor.
virtual Bool_t checkObservables(const RooArgSet *nset) const
Check if FUNC is valid for given normalization set.
void printMetaArgs(std::ostream &os) const
Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the p...
virtual void setCacheAndTrackHints(RooArgSet &)
Label OK'ed components of a RooRealSumPdf with cache-and-track.
Bool_t _extended
Definition: RooRealSumPdf.h:87
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise that all integrals can be handled internally.
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implement analytical integrations by deferring integration of component functions to integrators of c...
bool haveLastCoef() const
Definition: RooRealSumPdf.h:94
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
const Int_t n
Definition: legend1.C:16
@ InputArguments
Definition: RooGlobalFunc.h:68
Definition: first.py:1