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
22Implements 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
50#include "RooRealIntegral.h"
51#include "RooRealProxy.h"
52#include "RooRealVar.h"
53#include "RooMsgService.h"
54#include "RooNaNPacker.h"
55
56#include <TError.h>
57
58#include <algorithm>
59#include <memory>
60#include <stdexcept>
61
63
65
66////////////////////////////////////////////////////////////////////////////////
67/// Default constructor
68/// coverity[UNINIT_CTOR]
69
70RooRealSumPdf::RooRealSumPdf() : _normIntMgr(this, 10) {}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Constructor with name and title
74
75RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
76 RooAbsPdf(name,title),
77 _normIntMgr(this,10),
78 _funcList("!funcList","List of functions",this),
79 _coefList("!coefList","List of coefficients",this),
80 _extended(false)
81{
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Construct p.d.f consisting of \f$ \mathrm{coef}_1 * \mathrm{func}_1 + (1-\mathrm{coef}_1) * \mathrm{func}_2 \f$.
86/// The input coefficients and functions are allowed to be negative
87/// but the resulting sum is not, which is enforced at runtime.
88
89RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
90 RooAbsReal& func1, RooAbsReal& func2, RooAbsReal& coef1) :
91 RooRealSumPdf(name, title)
92{
93 // Special constructor with two functions and one coefficient
94
95 _funcList.add(func1) ;
96 _funcList.add(func2) ;
97 _coefList.add(coef1) ;
98
99}
100
101
102////////////////////////////////////////////////////////////////////////////////
103/// Constructor for a PDF from a list of functions and coefficients.
104/// It implements
105/// \f[
106/// \sum_i \mathrm{coef}_i \cdot \mathrm{func}_i,
107/// \f]
108/// if \f$ N_\mathrm{coef} = N_\mathrm{func} \f$. With `extended=true`, the coefficients can take any values. With `extended=false`,
109/// there is the danger of getting a degenerate minimisation problem because a PDF has to be normalised, which needs one degree
110/// of freedom less.
111///
112/// A plain (normalised) PDF can therefore be implemented with one less coefficient. RooFit then computes
113/// \f[
114/// \sum_i^{N-1} \mathrm{coef}_i \cdot \mathrm{func}_i + (1 - \sum_i \mathrm{coef}_i ) \cdot \mathrm{func}_N,
115/// \f]
116/// if \f$ N_\mathrm{coef} = N_\mathrm{func} - 1 \f$.
117///
118/// All coefficients and functions are allowed to be negative
119/// but the sum (*i.e.* the PDF) is not, which is enforced at runtime.
120///
121/// \param name Name of the PDF
122/// \param title Title (for plotting)
123/// \param inFuncList List of functions to sum
124/// \param inCoefList List of coefficients
125/// \param extended Interpret as extended PDF (requires equal number of functions and coefficients)
126
127RooRealSumPdf::RooRealSumPdf(const char *name, const char *title, const RooArgList &inFuncList,
128 const RooArgList &inCoefList, bool extended)
129 : RooRealSumPdf(name, title)
130{
131 setExtended(extended);
132
133 RooRealSumPdf::initializeFuncsAndCoefs(*this, inFuncList, inCoefList, _funcList, _coefList);
134}
135
136
138 const RooArgList& inFuncList, const RooArgList& inCoefList,
139 RooArgList& funcList, RooArgList& coefList)
140{
141 const std::string className = caller.ClassName();
142 const std::string constructorName = className + "::" + className;
143
144 if (!(inFuncList.size()==inCoefList.size()+1 || inFuncList.size()==inCoefList.size())) {
145 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName()
146 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << std::endl;
147 throw std::invalid_argument(className + ": 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 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") coefficient " << coef.GetName() << " is not of type RooAbsReal, ignored" << std::endl ;
157 continue ;
158 }
159 if (!dynamic_cast<const RooAbsReal*>(&func)) {
160 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") func " << func.GetName() << " is not of type RooAbsReal, ignored" << std::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 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") last func " << func.GetName() << " is not of type RooAbsReal, fatal error" << std::endl ;
171 throw std::invalid_argument(className + ": 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/// Collect the list of functions to be integrated from the cache
197
198const RooArgList &RooRealSumPdf::funcIntListFromCache(Int_t code, const char *rangeName) const
199{
200 return getCacheElem(*this, _normIntMgr, code, rangeName)->_funcIntList;
201}
202
204 Int_t code, const char *rangeName)
205{
206 return getCacheElem(caller, normIntMgr, code, rangeName)->_funcIntList;
207}
208
209////////////////////////////////////////////////////////////////////////////////
210
212{
214}
215
216////////////////////////////////////////////////////////////////////////////////
217
219 RooArgList const& funcList,
220 RooArgList const& coefList,
221 bool doFloor,
222 bool & hasWarnedBefore)
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 (!hasWarnedBefore) {
235 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
236 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
237 << sumCoeff << ". This means that the PDF is not properly normalised."
238 << " If the PDF was meant to be extended, provide as many coefficients as functions." << std::endl;
239 hasWarnedBefore = true;
240 }
241 // Signal that we are in an undefined region:
242 value = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
243 }
244
245 if (func->isSelectedComp()) {
246 value += func->getVal() * coefVal;
247 }
248
249 sumCoeff += coefVal;
250 }
251
252 // Introduce floor if so requested
253 return value < 0 && doFloor ? 0.0 : value;
254}
255
256////////////////////////////////////////////////////////////////////////////////
257/// Calculate the current value
258
260{
262}
263
264
266{
267 std::span<double> output = ctx.output();
268 std::size_t nEvents = output.size();
269
270 // Do running sum of coef/func pairs, calculate lastCoef.
271 std::fill(output.begin(), output.end(), 0.);
272
273 double sumCoeff = 0.;
274 for (unsigned int i = 0; i < _funcList.size(); ++i) {
275 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
276 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
277 const double coefVal = coef != nullptr ? ctx.at(coef)[0] : (1. - sumCoeff);
278
279 if (func->isSelectedComp()) {
280 auto funcValues = ctx.at(func);
281 if(funcValues.size() == 1) {
282 for (unsigned int j = 0; j < nEvents; ++j) {
283 output[j] += funcValues[0] * coefVal;
284 }
285 } else {
286 for (unsigned int j = 0; j < nEvents; ++j) {
287 output[j] += funcValues[j] * coefVal;
288 }
289 }
290 }
291
292 // Warn about degeneration of last coefficient
293 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
294 if (!_haveWarned) {
295 coutW(Eval) << "RooRealSumPdf::doEval(" << GetName()
296 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
297 << sumCoeff << ". This means that the PDF is not properly normalised. If the PDF was meant to be extended, provide as many coefficients as functions." << std::endl;
298 _haveWarned = true;
299 }
300 // Signal that we are in an undefined region by handing back one NaN.
301 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
302 }
303
304 sumCoeff += coefVal;
305 }
306
307 // Introduce floor if so requested
308 if (_doFloor || _doFloorGlobal) {
309 for (unsigned int j = 0; j < nEvents; ++j) {
310 output[j] = std::max(0., output[j]);
311 }
312 }
313}
314
316 RooArgList const& funcList, RooArgList const& coefList)
317{
318 bool ret(false) ;
319
320 for (unsigned int i=0; i < coefList.size(); ++i) {
321 const auto& coef = coefList[i];
322 const auto& func = funcList[i];
323
324 if (func.observableOverlaps(nset, coef)) {
325 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
326 << "): ERROR: coefficient " << coef.GetName()
327 << " and FUNC " << func.GetName() << " have one or more observables in common" << std::endl;
328 ret = true ;
329 }
330 if (coef.dependsOn(*nset)) {
331 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
332 << "): ERROR coefficient " << coef.GetName()
333 << " depends on one or more of the following observables" ; nset->Print("1") ;
334 ret = true ;
335 }
336 }
337
338 return ret ;
339}
340
341
342////////////////////////////////////////////////////////////////////////////////
343/// Check if FUNC is valid for given normalization set.
344/// Coefficient and FUNC must be non-overlapping, but func-coefficient
345/// pairs may overlap each other.
346///
347/// In the present implementation, coefficients may not be observables or derive
348/// from observables.
349
351{
352 return checkObservables(*this, nset, _funcList, _coefList);
353}
354
355
356////////////////////////////////////////////////////////////////////////////////
357/// Advertise that all integrals can be handled internally.
358
360 const RooArgSet* normSet2, const char* rangeName) const
361{
362 return getAnalyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, allVars, analVars, normSet2, rangeName);
363}
364
365
367 RooArgList const& funcList, RooArgList const& /*coefList*/,
368 RooArgSet& allVars, RooArgSet& analVars,
369 const RooArgSet* normSet2, const char* rangeName)
370{
371 // Handle trivial no-integration scenario
372 if (allVars.empty()) return 0 ;
373 if (caller.getForceNumInt()) return 0 ;
374
375 // Select subset of allVars that are actual dependents
376 analVars.add(allVars) ;
377 std::unique_ptr<RooArgSet> normSet;
378 if(normSet2) {
379 normSet = std::make_unique<RooArgSet>();
380 caller.getObservables(normSet2, *normSet);
381 }
382
383
384 // Check if this configuration was created before
385 Int_t sterileIdx(-1) ;
386 auto* cache = static_cast<CacheElem*>(normIntMgr.getObj(normSet.get(),&analVars,&sterileIdx,RooNameReg::ptr(rangeName)));
387 if (cache) {
388 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
389 return normIntMgr.lastIndex()+1 ;
390 }
391
392 // Create new cache element
393 cache = new CacheElem ;
394
395 // Make list of function projection and normalization integrals
396 for (const auto elm : funcList) {
397 const auto func = static_cast<RooAbsReal*>(elm);
398
399 std::unique_ptr<RooAbsReal> funcInt{func->createIntegral(analVars,rangeName)};
400 if(auto funcRealInt = dynamic_cast<RooRealIntegral*>(funcInt.get())) funcRealInt->setAllowComponentSelection(true);
401 cache->_funcIntList.addOwned(std::move(funcInt));
402 if (normSet && !normSet->empty()) {
403 cache->_funcNormList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(*normSet)});
404 }
405 }
406
407 // Store cache element
408 Int_t code = normIntMgr.setObj(normSet.get(),&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
409
410 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
411 return code+1 ;
412}
413
414////////////////////////////////////////////////////////////////////////////////
415/// Collect cache elements for given code to use in codegen
416
418RooRealSumPdf::getCacheElem(RooAbsReal const &caller, RooObjCacheManager &normIntMgr, Int_t code, const char *rangeName)
419{
420 // WVE needs adaptation for rangeName feature
421 auto *cache = static_cast<CacheElem *>(normIntMgr.getObjByIndex(code - 1));
422 if (cache == nullptr) { // revive the (sterilized) cache
423 // cout <<
424 // "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>")
425 // << ": reviving cache "<< endl;
426 RooArgSet vars;
427 caller.getParameters(nullptr, vars);
428 RooArgSet iset = normIntMgr.selectFromSet2(vars, code - 1);
429 RooArgSet nset = normIntMgr.selectFromSet1(vars, code - 1);
430 RooArgSet dummy;
431 Int_t code2 = caller.getAnalyticalIntegralWN(iset, dummy, &nset, rangeName);
432 R__ASSERT(code == code2); // must have revived the right (sterilized) slot...
433 cache = static_cast<CacheElem *>(normIntMgr.getObjByIndex(code - 1));
434 R__ASSERT(cache != nullptr);
435 }
436
437 return cache;
438}
439
440////////////////////////////////////////////////////////////////////////////////
441/// Implement analytical integrations by deferring integration of component
442/// functions to integrators of components.
443
444double RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
445{
446 return analyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, code, normSet2, rangeName, _haveWarned);
447}
448
450 RooArgList const& funcList, RooArgList const& coefList,
451 Int_t code, const RooArgSet* normSet2, const char* rangeName,
452 bool hasWarnedBefore)
453{
454 // Handle trivial passthrough scenario
455 if (code==0) return caller.getVal(normSet2) ;
456
457 // Retrieve cache element
458 const CacheElem *cache = getCacheElem(caller, normIntMgr, code, rangeName);
459
460 double value(0) ;
461
462 // N funcs, N-1 coefficients
463 double lastCoef(1) ;
464 auto funcIt = funcList.begin();
465 auto funcIntIt = cache->_funcIntList.begin();
466 for (const auto coefArg : coefList) {
467 assert(funcIt != funcList.end());
468 const auto coef = static_cast<const RooAbsReal*>(coefArg);
469 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
470 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
471
472 double coefVal = coef->getVal(normSet2) ;
473 if (coefVal) {
474 assert(func);
475 if (normSet2 ==nullptr || func->isSelectedComp()) {
476 assert(funcInt);
477 value += funcInt->getVal()*coefVal ;
478 }
479 lastCoef -= coef->getVal(normSet2) ;
480 }
481 }
482
483 const bool haveLastCoef = funcList.size() == coefList.size();
484
485 if (!haveLastCoef) {
486 // Add last func with correct coefficient
487 const auto func = static_cast<const RooAbsReal*>(*funcIt);
488 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
489 assert(func);
490
491 if (normSet2 ==nullptr || func->isSelectedComp()) {
492 assert(funcInt);
493 value += funcInt->getVal()*lastCoef ;
494 }
495
496 // Warn about coefficient degeneration
497 if (!hasWarnedBefore && (lastCoef<0 || lastCoef>1)) {
498 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
499 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
500 << 1-lastCoef << std::endl;
501 }
502 }
503
504 double normVal(1) ;
505 if (normSet2 && !normSet2->empty()) {
506 normVal = 0 ;
507
508 // N funcs, N-1 coefficients
509 auto funcNormIter = cache->_funcNormList.begin();
510 for (const auto coefAsArg : coefList) {
511 auto coef = static_cast<RooAbsReal*>(coefAsArg);
512 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
513
514 double coefVal = coef->getVal(normSet2);
515 if (coefVal) {
516 assert(funcNorm);
517 normVal += funcNorm->getVal()*coefVal ;
518 }
519 }
520
521 // Add last func with correct coefficient
522 if (!haveLastCoef) {
523 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
524 assert(funcNorm);
525
526 normVal += funcNorm->getVal()*lastCoef;
527 }
528 }
529
530 return value / normVal;
531}
532
533
534////////////////////////////////////////////////////////////////////////////////
535
537{
538 double n = getNorm(nset) ;
539 if (n<0) {
540 logEvalError("Expected number of events is negative") ;
541 }
542 return n ;
543}
544
545
546////////////////////////////////////////////////////////////////////////////////
547
548std::list<double>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
549{
550 return binBoundaries(_funcList, obs, xlo, xhi);
551}
552
553
554std::list<double> *
555RooRealSumPdf::binBoundaries(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
556{
557 std::unique_ptr<std::list<double>> sumBinB;
558 bool needClean(false);
559
560 // Loop over components pdf
561 for (auto *func : static_range_cast<RooAbsReal *>(funcList)) {
562
563 std::unique_ptr<std::list<double>> funcBinB{func->binBoundaries(obs, xlo, xhi)};
564
565 // Process hint
566 if (funcBinB) {
567 if (!sumBinB) {
568 // If this is the first hint, then just save it
569 sumBinB = std::move(funcBinB);
570 continue;
571 }
572
573 auto newSumBinB = std::make_unique<std::list<double>>(sumBinB->size() + funcBinB->size());
574
575 // Merge hints into temporary array
576 std::merge(funcBinB->begin(), funcBinB->end(), sumBinB->begin(), sumBinB->end(), newSumBinB->begin());
577
578 sumBinB = std::move(newSumBinB);
579 needClean = true;
580 }
581 }
582
583 // Remove consecutive duplicates
584 if (needClean) {
585 sumBinB->erase(std::unique(sumBinB->begin(), sumBinB->end()), sumBinB->end());
586 }
587
588 return sumBinB.release();
589}
590
591
592
593/// Check if all components that depend on `obs` are binned.
595{
596 return isBinnedDistribution(_funcList, obs);
597}
598
599
601{
602 for (auto* func : static_range_cast<RooAbsReal*>(funcList)) {
603
604 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
605 return false ;
606 }
607 }
608
609 return true ;
610}
611
612
613
614////////////////////////////////////////////////////////////////////////////////
615
616std::list<double>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
617{
618 return plotSamplingHint(_funcList, obs, xlo, xhi);
619}
620
621std::list<double> *
622RooRealSumPdf::plotSamplingHint(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
623{
624 std::unique_ptr<std::list<double>> sumHint;
625 bool needClean(false);
626
627 // Loop over components pdf
628 for (const auto elm : funcList) {
629 auto func = static_cast<RooAbsReal *>(elm);
630
631 std::unique_ptr<std::list<double>> funcHint{func->plotSamplingHint(obs, xlo, xhi)};
632
633 // Process hint
634 if (funcHint) {
635 if (!sumHint) {
636
637 // If this is the first hint, then just save it
638 sumHint = std::move(funcHint);
639 continue;
640 }
641
642 auto newSumHint = std::make_unique<std::list<double>>(sumHint->size() + funcHint->size());
643
644 // the lists must be sorted before merging them
645 funcHint->sort();
646 sumHint->sort();
647 // Merge hints into temporary array
648 std::merge(funcHint->begin(), funcHint->end(), sumHint->begin(), sumHint->end(), newSumHint->begin());
649
650 sumHint = std::move(newSumHint);
651 needClean = true;
652 }
653 }
654
655 // Remove consecutive duplicates
656 if (needClean) {
657 sumHint->erase(std::unique(sumHint->begin(), sumHint->end()), sumHint->end());
658 }
659
660 return sumHint.release();
661}
662
663
664
665
666////////////////////////////////////////////////////////////////////////////////
667/// Label OK'ed components of a RooRealSumPdf with cache-and-track
668
670{
671 setCacheAndTrackHints(_funcList, trackNodes);
672}
673
674
676{
677 for (const auto sarg : funcList) {
678 if (sarg->canNodeBeCached()==Always) {
679 trackNodes.add(*sarg) ;
680 }
681 }
682}
683
684
685////////////////////////////////////////////////////////////////////////////////
686/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
687/// product operator construction
688
689void RooRealSumPdf::printMetaArgs(std::ostream& os) const
690{
692}
693
694
695void RooRealSumPdf::printMetaArgs(RooArgList const& funcList, RooArgList const& coefList, std::ostream& os)
696{
697
698 bool first(true) ;
699
700 if (!coefList.empty()) {
701 auto funcIter = funcList.begin();
702
703 for (const auto coef : coefList) {
704 if (!first) {
705 os << " + " ;
706 } else {
707 first = false ;
708 }
709 const auto func = *(funcIter++);
710 os << coef->GetName() << " * " << func->GetName();
711 }
712
713 if (funcIter != funcList.end()) {
714 os << " + [%] * " << (*funcIter)->GetName() ;
715 }
716 } else {
717
718 for (const auto func : funcList) {
719 if (!first) {
720 os << " + " ;
721 } else {
722 first = false ;
723 }
724 os << func->GetName() ;
725 }
726 }
727
728 os << " " ;
729}
730
731std::unique_ptr<RooAbsArg>
733{
734 if (normSet.empty() || selfNormalized()) {
735 return RooAbsPdf::compileForNormSet({}, ctx);
736 }
737 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
738
739 if (ctx.likelihoodMode() && pdfClone->getAttribute("BinnedLikelihood")) {
740
741 // This has to be done before compiling the servers, such that the
742 // RooBinWidthFunctions know to disable themselves.
743 ctx.setBinnedLikelihoodMode(true);
744
745 ctx.markAsCompiled(*pdfClone);
746 ctx.compileServers(*pdfClone, {});
747
748 pdfClone->setAttribute("BinnedLikelihoodActive");
749 // If this is a binned likelihood, we're flagging it in the context.
750 // Then, the RooBinWidthFunctions know that they should not put
751 // themselves in the computation graph. Like this, the pdf values can
752 // directly be interpreted as yields, without multiplying them with the
753 // bin widths again in the NLL. However, the NLL class has to be careful
754 // to only skip the bin with multiplication when there actually were
755 // RooBinWidthFunctions! This is not the case for old workspace before
756 // ROOT 6.26. Therefore, we use the "BinnedLikelihoodActiveYields"
757 // attribute to let the NLL know what it should do.
758 if (ctx.binWidthFuncFlag()) {
759 pdfClone->setAttribute("BinnedLikelihoodActiveYields");
760 }
761 return pdfClone;
762 }
763
764 ctx.compileServers(*pdfClone, {});
765
766 RooArgSet depList;
767 pdfClone->getObservables(&normSet, depList);
768
769 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, depList);
770
771 // The direct servers are this pdf and the normalization integral, which
772 // don't need to be compiled further.
773 for (RooAbsArg *server : newArg->servers()) {
774 ctx.markAsCompiled(*server);
775 }
776 ctx.markAsCompiled(*newArg);
777 newArg->addOwnedComponents(std::move(pdfClone));
778 return newArg;
779}
780
781std::unique_ptr<RooAbsReal> RooRealSumPdf::createExpectedEventsFunc(const RooArgSet *nset) const
782{
783 if (nset == nullptr)
784 return nullptr;
785 return std::unique_ptr<RooAbsReal>{createIntegral(*nset, *getIntegratorConfig(), normRange())};
786}
#define oocoutW(o, a)
#define coutW(a)
#define oocoutE(o, a)
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooFit::OwningPtr< 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...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:89
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
const_iterator end() const
Storage_t::size_type size() const
const_iterator begin() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:195
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
@ CanBeExtended
Definition RooAbsPdf.h:212
@ CanNotBeExtended
Definition RooAbsPdf.h:212
const char * normRange() const
Definition RooAbsPdf.h:250
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
bool getForceNumInt() const
Definition RooAbsReal.h:176
void logEvalError(const char *message, const char *serverValueString=nullptr) const
Log evaluation error message.
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
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...
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
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:24
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
RooArgSet selectFromSet1(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 1 with a given index and an i...
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 2 with a given index and an i...
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
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...
void markAsCompiled(RooAbsArg &arg) const
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
Implements a PDF constructed from a sum of functions:
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components of a RooRealSumPdf with cache-and-track.
static bool _doFloorGlobal
Global flag for introducing floor at zero in pdf.
bool _doFloor
Introduce floor at zero in pdf.
std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const override
Returns an object that represents the expected number of events for a given normalization set,...
bool checkObservables(const RooArgSet *nset) const override
Check if FUNC is valid for given normalization set.
void setExtended(bool extended)
const RooArgList & funcList() const
RooObjCacheManager _normIntMgr
! The integration cache manager
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
std::list< double > * plotSamplingHint(RooAbsRealLValue &, double, double) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
double evaluate() const override
Calculate the current value.
ExtendMode extendMode() const override
Returns ability of PDF to provide extended likelihood terms.
const RooArgList & funcIntListFromCache(Int_t code, const char *rangeName=nullptr) const
Collect the list of functions to be integrated from the cache.
RooListProxy _coefList
List of coefficients.
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
RooListProxy _funcList
List of component FUNCs.
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Retrieve bin boundaries if this distribution is binned in obs.
RooRealSumPdf()
Default constructor coverity[UNINIT_CTOR].
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the p...
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
double expectedEvents(const RooArgSet *nset) const override
Return expected number of events for extended likelihood calculation, which is the sum of all coeffic...
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Implement analytical integrations by deferring integration of component functions to integrators of c...
static const CacheElem * getCacheElem(RooAbsReal const &caller, RooObjCacheManager &normIntMgr, Int_t code, const char *rangeName)
Collect cache elements for given code to use in codegen.
static void initializeFuncsAndCoefs(RooAbsReal const &caller, const RooArgList &inFuncList, const RooArgList &inCoefList, RooArgList &funcList, RooArgList &coefList)
const RooArgList & coefList() const
bool _extended
Allow use as extended p.d.f.
bool isBinnedDistribution(const RooArgSet &obs) const override
Check if all components that depend on obs are binned.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
const Int_t n
Definition legend1.C:16
__roodevice__ static __roohost__ double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output()