Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "RooRealSumPdf.h"
48
49#include "RooRealIntegral.h"
50#include "RooRealProxy.h"
51#include "RooRealVar.h"
52#include "RooMsgService.h"
53#include "RooNaNPacker.h"
54#include "RunContext.h"
55
56#include <TError.h>
57
58#include <algorithm>
59#include <memory>
60#include <stdexcept>
61
62using namespace std;
63
65
67
68////////////////////////////////////////////////////////////////////////////////
69/// Default constructor
70/// coverity[UNINIT_CTOR]
71
72RooRealSumPdf::RooRealSumPdf() : _normIntMgr(this,10)
73{
74
75}
76
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// Constructor with name and title
81
82RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
83 RooAbsPdf(name,title),
84 _normIntMgr(this,10),
85 _funcList("!funcList","List of functions",this),
86 _coefList("!coefList","List of coefficients",this),
87 _extended(false),
88 _doFloor(false)
89{
90
91}
92
93
94
95////////////////////////////////////////////////////////////////////////////////
96/// Construct p.d.f consisting of \f$ \mathrm{coef}_1 * \mathrm{func}_1 + (1-\mathrm{coef}_1) * \mathrm{func}_2 \f$.
97/// The input coefficients and functions are allowed to be negative
98/// but the resulting sum is not, which is enforced at runtime.
99
100RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
101 RooAbsReal& func1, RooAbsReal& func2, RooAbsReal& coef1) :
102 RooRealSumPdf(name, title)
103{
104 // Special constructor with two functions and one coefficient
105
106 _funcList.add(func1) ;
107 _funcList.add(func2) ;
108 _coefList.add(coef1) ;
109
110}
111
112
113////////////////////////////////////////////////////////////////////////////////
114/// Constructor for a PDF from a list of functions and coefficients.
115/// It implements
116/// \f[
117/// \sum_i \mathrm{coef}_i \cdot \mathrm{func}_i,
118/// \f]
119/// if \f$ N_\mathrm{coef} = N_\mathrm{func} \f$. With `extended=true`, the coefficients can take any values. With `extended=false`,
120/// there is the danger of getting a degenerate minimisation problem because a PDF has to be normalised, which needs one degree
121/// of freedom less.
122///
123/// A plain (normalised) PDF can therefore be implemented with one less coefficient. RooFit then computes
124/// \f[
125/// \sum_i^{N-1} \mathrm{coef}_i \cdot \mathrm{func}_i + (1 - \sum_i \mathrm{coef}_i ) \cdot \mathrm{func}_N,
126/// \f]
127/// if \f$ N_\mathrm{coef} = N_\mathrm{func} - 1 \f$.
128///
129/// All coefficients and functions are allowed to be negative
130/// but the sum (*i.e.* the PDF) is not, which is enforced at runtime.
131///
132/// \param name Name of the PDF
133/// \param title Title (for plotting)
134/// \param inFuncList List of functions to sum
135/// \param inCoefList List of coefficients
136/// \param extended Interpret as extended PDF (requires equal number of functions and coefficients)
137
138RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
139 const RooArgList& inFuncList, const RooArgList& inCoefList, Bool_t extended) :
140 RooRealSumPdf(name, title)
141{
142 _extended = extended;
143
144 if (!(inFuncList.getSize()==inCoefList.getSize()+1 || inFuncList.getSize()==inCoefList.getSize())) {
145 coutE(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName()
146 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << endl ;
147 throw std::invalid_argument("RooRealSumPdf: Number of PDFs and coefficients is inconsistent.");
148 }
149
150 // Constructor with N functions and N or N-1 coefs
151 for (unsigned int i = 0; i < inCoefList.size(); ++i) {
152 const auto& func = inFuncList[i];
153 const auto& coef = inCoefList[i];
154
155 if (!dynamic_cast<const RooAbsReal*>(&coef)) {
156 coutW(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") coefficient " << coef.GetName() << " is not of type RooAbsReal, ignored" << endl ;
157 continue ;
158 }
159 if (!dynamic_cast<const RooAbsReal*>(&func)) {
160 coutW(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") func " << func.GetName() << " is not of type RooAbsReal, ignored" << endl ;
161 continue ;
162 }
163 _funcList.add(func) ;
164 _coefList.add(coef) ;
165 }
166
167 if (inFuncList.size() == inCoefList.size() + 1) {
168 const auto& func = inFuncList[inFuncList.size()-1];
169 if (!dynamic_cast<const RooAbsReal*>(&func)) {
170 coutE(InputArguments) << "RooRealSumPdf::RooRealSumPdf(" << GetName() << ") last func " << func.GetName() << " is not of type RooAbsReal, fatal error" << endl ;
171 throw std::invalid_argument("RooRealSumPdf: Function passed as is not of type RooAbsReal.");
172 }
173 _funcList.add(func);
174 }
175
176}
177
178
179
180
181////////////////////////////////////////////////////////////////////////////////
182/// Copy constructor
183
185 RooAbsPdf(other,name),
186 _normIntMgr(other._normIntMgr,this),
187 _funcList("!funcList",this,other._funcList),
188 _coefList("!coefList",this,other._coefList),
189 _extended(other._extended),
190 _doFloor(other._doFloor)
191{
192
193}
194
195
196
197////////////////////////////////////////////////////////////////////////////////
198/// Destructor
199
201{
202
203}
204
205
206
207
208
209////////////////////////////////////////////////////////////////////////////////
210
212{
214}
215
216
217
218
219////////////////////////////////////////////////////////////////////////////////
220/// Calculate the current value
221
223{
224 // Do running sum of coef/func pairs, calculate lastCoef.
225 double value = 0;
226 double sumCoeff = 0.;
227 for (unsigned int i = 0; i < _funcList.size(); ++i) {
228 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
229 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
230 const double coefVal = coef != nullptr ? coef->getVal() : (1. - sumCoeff);
231
232 // Warn about degeneration of last coefficient
233 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
234 if (!_haveWarned) {
235 coutW(Eval) << "RooRealSumPdf::evaluate(" << GetName()
236 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
237 << sumCoeff << ". This means that the PDF is not properly normalised. If the PDF was meant to be extended, provide as many coefficients as functions." << endl ;
238 _haveWarned = true;
239 }
240 // Signal that we are in an undefined region:
241 value = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
242 }
243
244 if (func->isSelectedComp()) {
245 value += func->getVal() * coefVal;
246 }
247
248 sumCoeff += coefVal;
249 }
250
251 // Introduce floor if so requested
252 if (value<0 && (_doFloor || _doFloorGlobal)) {
253 value = 0 ;
254 }
255
256 return value ;
257}
258
259
260void RooRealSumPdf::computeBatch(cudaStream_t* /*stream*/, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const {
261 // Do running sum of coef/func pairs, calculate lastCoef.
262 for (unsigned int j = 0; j < nEvents; ++j) {
263 output[j] = 0.0;
264 }
265
266 double sumCoeff = 0.;
267 for (unsigned int i = 0; i < _funcList.size(); ++i) {
268 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
269 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
270 const double coefVal = coef != nullptr ? dataMap.at(coef)[0] : (1. - sumCoeff);
271
272 if (func->isSelectedComp()) {
273 auto funcValues = dataMap.at(func);
274 if(funcValues.size() == 1) {
275 for (unsigned int j = 0; j < nEvents; ++j) {
276 output[j] += funcValues[0] * coefVal;
277 }
278 } else {
279 for (unsigned int j = 0; j < nEvents; ++j) {
280 output[j] += funcValues[j] * coefVal;
281 }
282 }
283 }
284
285 // Warn about degeneration of last coefficient
286 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
287 if (!_haveWarned) {
288 coutW(Eval) << "RooRealSumPdf::evaluateSpan(" << GetName()
289 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
290 << sumCoeff << ". This means that the PDF is not properly normalised. If the PDF was meant to be extended, provide as many coefficients as functions." << endl ;
291 _haveWarned = true;
292 }
293 // Signal that we are in an undefined region by handing back one NaN.
294 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
295 }
296
297 sumCoeff += coefVal;
298 }
299
300 // Introduce floor if so requested
301 if (_doFloor || _doFloorGlobal) {
302 for (unsigned int j = 0; j < nEvents; ++j) {
303 output[j] += std::max(0., output[j]);
304 }
305 }
306}
307
308
309////////////////////////////////////////////////////////////////////////////////
310/// Check if FUNC is valid for given normalization set.
311/// Coefficient and FUNC must be non-overlapping, but func-coefficient
312/// pairs may overlap each other.
313///
314/// In the present implementation, coefficients may not be observables or derive
315/// from observables.
316
318{
319 Bool_t ret(kFALSE) ;
320
321 for (unsigned int i=0; i < _coefList.size(); ++i) {
322 const auto& coef = _coefList[i];
323 const auto& func = _funcList[i];
324
325 if (func.observableOverlaps(nset, coef)) {
326 coutE(InputArguments) << "RooRealSumPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef.GetName()
327 << " and FUNC " << func.GetName() << " have one or more observables in common" << endl ;
328 ret = kTRUE ;
329 }
330 if (coef.dependsOn(*nset)) {
331 coutE(InputArguments) << "RooRealPdf::checkObservables(" << GetName() << "): ERROR coefficient " << coef.GetName()
332 << " depends on one or more of the following observables" ; nset->Print("1") ;
333 ret = kTRUE ;
334 }
335 }
336
337 return ret ;
338}
339
340
341
342
343////////////////////////////////////////////////////////////////////////////////
344/// Advertise that all integrals can be handled internally.
345
347 const RooArgSet* normSet2, const char* rangeName) const
348{
349 // Handle trivial no-integration scenario
350 if (allVars.getSize()==0) return 0 ;
351 if (_forceNumInt) return 0 ;
352
353 // Select subset of allVars that are actual dependents
354 analVars.add(allVars) ;
355 RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
356
357
358 // Check if this configuration was created before
359 Int_t sterileIdx(-1) ;
360 CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx,RooNameReg::ptr(rangeName)) ;
361 if (cache) {
362 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
363 return _normIntMgr.lastIndex()+1 ;
364 }
365
366 // Create new cache element
367 cache = new CacheElem ;
368
369 // Make list of function projection and normalization integrals
370 for (const auto elm : _funcList) {
371 const auto func = static_cast<RooAbsReal*>(elm);
372
373 RooAbsReal* funcInt = func->createIntegral(analVars,rangeName) ;
374 if(funcInt->InheritsFrom(RooRealIntegral::Class())) ((RooRealIntegral*)funcInt)->setAllowComponentSelection(true);
375 cache->_funcIntList.addOwned(*funcInt) ;
376 if (normSet && normSet->getSize()>0) {
377 RooAbsReal* funcNorm = func->createIntegral(*normSet) ;
378 cache->_funcNormList.addOwned(*funcNorm) ;
379 }
380 }
381
382 // Store cache element
383 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
384
385 if (normSet) {
386 delete normSet ;
387 }
388
389 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
390 return code+1 ;
391}
392
393
394
395
396////////////////////////////////////////////////////////////////////////////////
397/// Implement analytical integrations by deferring integration of component
398/// functions to integrators of components.
399
400Double_t RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
401{
402 // Handle trivial passthrough scenario
403 if (code==0) return getVal(normSet2) ;
404
405
406 // WVE needs adaptation for rangeName feature
407 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
408 if (cache==0) { // revive the (sterilized) cache
409 //cout << "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << ": reviving cache "<< endl;
410 std::unique_ptr<RooArgSet> vars( getParameters(RooArgSet()) );
411 RooArgSet iset = _normIntMgr.selectFromSet2(*vars, code-1);
412 RooArgSet nset = _normIntMgr.selectFromSet1(*vars, code-1);
413 RooArgSet dummy;
414 Int_t code2 = getAnalyticalIntegralWN(iset,dummy,&nset,rangeName);
415 R__ASSERT(code==code2); // must have revived the right (sterilized) slot...
416 cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
417 R__ASSERT(cache!=0);
418 }
419
420 Double_t value(0) ;
421
422 // N funcs, N-1 coefficients
423 Double_t lastCoef(1) ;
424 auto funcIt = _funcList.begin();
425 auto funcIntIt = cache->_funcIntList.begin();
426 for (const auto coefArg : _coefList) {
427 assert(funcIt != _funcList.end());
428 const auto coef = static_cast<const RooAbsReal*>(coefArg);
429 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
430 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
431
432 Double_t coefVal = coef->getVal(normSet2) ;
433 if (coefVal) {
434 assert(func);
435 if (normSet2 ==0 || func->isSelectedComp()) {
436 assert(funcInt);
437 value += funcInt->getVal()*coefVal ;
438 }
439 lastCoef -= coef->getVal(normSet2) ;
440 }
441 }
442
443 if (!haveLastCoef()) {
444 // Add last func with correct coefficient
445 const auto func = static_cast<const RooAbsReal*>(*funcIt);
446 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
447 assert(func);
448
449 if (normSet2 ==0 || func->isSelectedComp()) {
450 assert(funcInt);
451 value += funcInt->getVal()*lastCoef ;
452 }
453
454 // Warn about coefficient degeneration
455 if (!_haveWarned && (lastCoef<0 || lastCoef>1)) {
456 coutW(Eval) << "RooRealSumPdf::evaluate(" << GetName()
457 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
458 << 1-lastCoef << endl ;
459 }
460 }
461
462 Double_t normVal(1) ;
463 if (normSet2 && normSet2->getSize()>0) {
464 normVal = 0 ;
465
466 // N funcs, N-1 coefficients
467 auto funcNormIter = cache->_funcNormList.begin();
468 for (const auto coefAsArg : _coefList) {
469 auto coef = static_cast<RooAbsReal*>(coefAsArg);
470 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
471
472 Double_t coefVal = coef->getVal(normSet2);
473 if (coefVal) {
474 assert(funcNorm);
475 normVal += funcNorm->getVal()*coefVal ;
476 }
477 }
478
479 // Add last func with correct coefficient
480 if (!haveLastCoef()) {
481 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
482 assert(funcNorm);
483
484 normVal += funcNorm->getVal()*lastCoef;
485 }
486 }
487
488 return value / normVal;
489}
490
491
492////////////////////////////////////////////////////////////////////////////////
493
495{
496 Double_t n = getNorm(nset) ;
497 if (n<0) {
498 logEvalError("Expected number of events is negative") ;
499 }
500 return n ;
501}
502
503
504////////////////////////////////////////////////////////////////////////////////
505
506std::list<Double_t>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
507{
508 list<Double_t>* sumBinB = 0 ;
509 Bool_t needClean(kFALSE) ;
510
511 // Loop over components pdf
512 for (const auto elm : _funcList) {
513 auto func = static_cast<RooAbsReal*>(elm);
514
515 list<Double_t>* funcBinB = func->binBoundaries(obs,xlo,xhi) ;
516
517 // Process hint
518 if (funcBinB) {
519 if (!sumBinB) {
520 // If this is the first hint, then just save it
521 sumBinB = funcBinB ;
522 } else {
523
524 list<Double_t>* newSumBinB = new list<Double_t>(sumBinB->size()+funcBinB->size()) ;
525
526 // Merge hints into temporary array
527 merge(funcBinB->begin(),funcBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
528
529 // Copy merged array without duplicates to new sumBinBArrau
530 delete sumBinB ;
531 delete funcBinB ;
532 sumBinB = newSumBinB ;
533 needClean = kTRUE ;
534 }
535 }
536 }
537
538 // Remove consecutive duplicates
539 if (needClean) {
540 list<Double_t>::iterator new_end = unique(sumBinB->begin(),sumBinB->end()) ;
541 sumBinB->erase(new_end,sumBinB->end()) ;
542 }
543
544 return sumBinB ;
545}
546
547
548
549/// Check if all components that depend on `obs` are binned.
551{
552 for (const auto elm : _funcList) {
553 auto func = static_cast<RooAbsReal*>(elm);
554
555 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
556 return kFALSE ;
557 }
558 }
559
560 return kTRUE ;
561}
562
563
564
565
566
567////////////////////////////////////////////////////////////////////////////////
568
569std::list<Double_t>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
570{
571 list<Double_t>* sumHint = 0 ;
572 Bool_t needClean(kFALSE) ;
573
574 // Loop over components pdf
575 for (const auto elm : _funcList) {
576 auto func = static_cast<RooAbsReal*>(elm);
577
578 list<Double_t>* funcHint = func->plotSamplingHint(obs,xlo,xhi) ;
579
580 // Process hint
581 if (funcHint) {
582 if (!sumHint) {
583
584 // If this is the first hint, then just save it
585 sumHint = funcHint ;
586
587 } else {
588
589 list<Double_t>* newSumHint = new list<Double_t>(sumHint->size()+funcHint->size()) ;
590
591 // Merge hints into temporary array
592 merge(funcHint->begin(),funcHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
593
594 // Copy merged array without duplicates to new sumHintArrau
595 delete sumHint ;
596 sumHint = newSumHint ;
597 needClean = kTRUE ;
598 }
599 }
600 }
601
602 // Remove consecutive duplicates
603 if (needClean) {
604 list<Double_t>::iterator new_end = unique(sumHint->begin(),sumHint->end()) ;
605 sumHint->erase(new_end,sumHint->end()) ;
606 }
607
608 return sumHint ;
609}
610
611
612
613
614////////////////////////////////////////////////////////////////////////////////
615/// Label OK'ed components of a RooRealSumPdf with cache-and-track
616
618{
619 for (const auto sarg : _funcList) {
620 if (sarg->canNodeBeCached()==Always) {
621 trackNodes.add(*sarg) ;
622 //cout << "tracking node RealSumPdf component " << sarg->IsA()->GetName() << "::" << sarg->GetName() << endl ;
623 }
624 }
625}
626
627
628
629////////////////////////////////////////////////////////////////////////////////
630/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
631/// product operator construction
632
633void RooRealSumPdf::printMetaArgs(ostream& os) const
634{
635
637
638 if (_coefList.getSize()!=0) {
639 auto funcIter = _funcList.begin();
640
641 for (const auto coef : _coefList) {
642 if (!first) {
643 os << " + " ;
644 } else {
645 first = kFALSE ;
646 }
647 const auto func = *(funcIter++);
648 os << coef->GetName() << " * " << func->GetName();
649 }
650
651 if (funcIter != _funcList.end()) {
652 os << " + [%] * " << (*funcIter)->GetName() ;
653 }
654 } else {
655
656 for (const auto func : _funcList) {
657 if (!first) {
658 os << " + " ;
659 } else {
660 first = kFALSE ;
661 }
662 os << func->GetName() ;
663 }
664 }
665
666 os << " " ;
667}
#define coutW(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:101
bool Bool_t
Definition RtypesCore.h:63
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
#define R__ASSERT(e)
Definition TError.h:118
char name[80]
Definition TGX11.cxx:110
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
friend class RooArgSet
Definition RooAbsArg.h:642
RooArgSet * getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
const_iterator end() const
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add an argument and transfer the ownership to the collection.
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.
const char * GetName() const
Returns name of object.
Double_t getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:239
@ CanBeExtended
Definition RooAbsPdf.h:256
@ CanNotBeExtended
Definition RooAbsPdf.h:256
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:64
Bool_t _forceNumInt
Definition RooAbsReal.h:487
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 ...
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.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
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...
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:22
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)
RooArgSet selectFromSet1(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 1.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 2.
Int_t lastIndex() const
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
RooObjCacheManager _normIntMgr
Double_t evaluate() const
Calculate the current value.
RooListProxy _coefList
RooListProxy _funcList
The integration cache manager.
RooRealSumPdf()
Default constructor coverity[UNINIT_CTOR].
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events for extended likelihood calculation, which is the sum of all coeffic...
Bool_t isBinnedDistribution(const RooArgSet &obs) const
Check if all components that depend on obs are binned.
static Bool_t _doFloorGlobal
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &, Double_t, Double_t) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
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.
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Retrieve bin boundaries if this distribution is binned in obs.
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
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
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:515
const Int_t n
Definition legend1.C:16
Definition first.py:1
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output(int code)
Definition gifencode.c:226