Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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,
91 RooRealSumPdf(name, title)
92{
93 // Special constructor with two functions and one coefficient
94
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
135
136
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
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
201
202
204 RooArgList const& funcList,
205 RooArgList const& coefList,
206 bool doFloor,
207 bool & hasWarnedBefore)
208{
209 // Do running sum of coef/func pairs, calculate lastCoef.
210 double value = 0;
211 double sumCoeff = 0.;
212 for (unsigned int i = 0; i < funcList.size(); ++i) {
213 const auto func = static_cast<RooAbsReal*>(&funcList[i]);
214 const auto coef = static_cast<RooAbsReal*>(i < coefList.size() ? &coefList[i] : nullptr);
215 const double coefVal = coef != nullptr ? coef->getVal() : (1. - sumCoeff);
216
217 // Warn about degeneration of last coefficient
218 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
219 if (!hasWarnedBefore) {
220 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
221 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
222 << sumCoeff << ". This means that the PDF is not properly normalised."
223 << " If the PDF was meant to be extended, provide as many coefficients as functions." << std::endl;
224 hasWarnedBefore = true;
225 }
226 // Signal that we are in an undefined region:
227 value = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
228 }
229
230 if (func->isSelectedComp()) {
231 value += func->getVal() * coefVal;
232 }
233
234 sumCoeff += coefVal;
235 }
236
237 // Introduce floor if so requested
238 return value < 0 && doFloor ? 0.0 : value;
239}
240
241////////////////////////////////////////////////////////////////////////////////
242/// Calculate the current value
243
245{
247}
248
249
251{
252 std::span<double> output = ctx.output();
253 std::size_t nEvents = output.size();
254
255 // Do running sum of coef/func pairs, calculate lastCoef.
256 std::fill(output.begin(), output.end(), 0.);
257
258 double sumCoeff = 0.;
259 for (unsigned int i = 0; i < _funcList.size(); ++i) {
260 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
261 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
262 const double coefVal = coef != nullptr ? ctx.at(coef)[0] : (1. - sumCoeff);
263
264 if (func->isSelectedComp()) {
265 auto funcValues = ctx.at(func);
266 if(funcValues.size() == 1) {
267 for (unsigned int j = 0; j < nEvents; ++j) {
268 output[j] += funcValues[0] * coefVal;
269 }
270 } else {
271 for (unsigned int j = 0; j < nEvents; ++j) {
272 output[j] += funcValues[j] * coefVal;
273 }
274 }
275 }
276
277 // Warn about degeneration of last coefficient
278 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
279 if (!_haveWarned) {
280 coutW(Eval) << "RooRealSumPdf::doEval(" << GetName()
281 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
282 << 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;
283 _haveWarned = true;
284 }
285 // Signal that we are in an undefined region by handing back one NaN.
286 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
287 }
288
289 sumCoeff += coefVal;
290 }
291
292 // Introduce floor if so requested
293 if (_doFloor || _doFloorGlobal) {
294 for (unsigned int j = 0; j < nEvents; ++j) {
295 output[j] = std::max(0., output[j]);
296 }
297 }
298}
299
301 RooArgList const &funcList, RooArgList const &coefList, bool normalize)
302{
303 bool noLastCoeff = funcList.size() != coefList.size();
304
305 std::string const &funcName = ctx.buildArg(funcList);
306 std::string const &coeffName = ctx.buildArg(coefList);
307 std::string const &coeffSize = std::to_string(coefList.size());
308
309 std::string sum = ctx.getTmpVarName();
310 std::string coeffSum = ctx.getTmpVarName();
311 ctx.addToCodeBody(klass, "double " + sum + " = 0;\ndouble " + coeffSum + "= 0;\n");
312
313 std::string iterator = "i_" + ctx.getTmpVarName();
314 std::string subscriptExpr = "[" + iterator + "]";
315
316 std::string code = "for(int " + iterator + " = 0; " + iterator + " < " + coeffSize + "; " + iterator + "++) {\n" +
317 sum + " += " + funcName + subscriptExpr + " * " + coeffName + subscriptExpr + ";\n";
318 code += coeffSum + " += " + coeffName + subscriptExpr + ";\n";
319 code += "}\n";
320
321 if (noLastCoeff) {
322 code += sum + " += " + funcName + "[" + coeffSize + "]" + " * (1 - " + coeffSum + ");\n";
323 } else if (normalize) {
324 code += sum + " /= " + coeffSum + ";\n";
325 }
326 ctx.addToCodeBody(klass, code);
327
328 return sum;
329}
330
335
337 RooArgList const& funcList, RooArgList const& coefList)
338{
339 bool ret(false) ;
340
341 for (unsigned int i=0; i < coefList.size(); ++i) {
342 const auto& coef = coefList[i];
343 const auto& func = funcList[i];
344
345 if (func.observableOverlaps(nset, coef)) {
346 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
347 << "): ERROR: coefficient " << coef.GetName()
348 << " and FUNC " << func.GetName() << " have one or more observables in common" << std::endl;
349 ret = true ;
350 }
351 if (coef.dependsOn(*nset)) {
352 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
353 << "): ERROR coefficient " << coef.GetName()
354 << " depends on one or more of the following observables" ; nset->Print("1") ;
355 ret = true ;
356 }
357 }
358
359 return ret ;
360}
361
362
363////////////////////////////////////////////////////////////////////////////////
364/// Check if FUNC is valid for given normalization set.
365/// Coefficient and FUNC must be non-overlapping, but func-coefficient
366/// pairs may overlap each other.
367///
368/// In the present implementation, coefficients may not be observables or derive
369/// from observables.
370
372{
373 return checkObservables(*this, nset, _funcList, _coefList);
374}
375
376
377////////////////////////////////////////////////////////////////////////////////
378/// Advertise that all integrals can be handled internally.
379
381 const RooArgSet* normSet2, const char* rangeName) const
382{
383 return getAnalyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, allVars, analVars, normSet2, rangeName);
384}
385
386
388 RooArgList const& funcList, RooArgList const& /*coefList*/,
389 RooArgSet& allVars, RooArgSet& analVars,
390 const RooArgSet* normSet2, const char* rangeName)
391{
392 // Handle trivial no-integration scenario
393 if (allVars.empty()) return 0 ;
394 if (caller.getForceNumInt()) return 0 ;
395
396 // Select subset of allVars that are actual dependents
397 analVars.add(allVars) ;
398 std::unique_ptr<RooArgSet> normSet;
399 if(normSet2) {
400 normSet = std::make_unique<RooArgSet>();
401 caller.getObservables(normSet2, *normSet);
402 }
403
404
405 // Check if this configuration was created before
406 Int_t sterileIdx(-1) ;
407 auto* cache = static_cast<CacheElem*>(normIntMgr.getObj(normSet.get(),&analVars,&sterileIdx,RooNameReg::ptr(rangeName)));
408 if (cache) {
409 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
410 return normIntMgr.lastIndex()+1 ;
411 }
412
413 // Create new cache element
414 cache = new CacheElem ;
415
416 // Make list of function projection and normalization integrals
417 for (const auto elm : funcList) {
418 const auto func = static_cast<RooAbsReal*>(elm);
419
420 std::unique_ptr<RooAbsReal> funcInt{func->createIntegral(analVars,rangeName)};
421 if(auto funcRealInt = dynamic_cast<RooRealIntegral*>(funcInt.get())) funcRealInt->setAllowComponentSelection(true);
422 cache->_funcIntList.addOwned(std::move(funcInt));
423 if (normSet && !normSet->empty()) {
424 cache->_funcNormList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(*normSet)});
425 }
426 }
427
428 // Store cache element
429 Int_t code = normIntMgr.setObj(normSet.get(),&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
430
431 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
432 return code+1 ;
433}
434
435
436
437
438////////////////////////////////////////////////////////////////////////////////
439/// Implement analytical integrations by deferring integration of component
440/// functions to integrators of components.
441
442double RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
443{
444 return analyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, code, normSet2, rangeName, _haveWarned);
445}
446
448 RooArgList const& funcList, RooArgList const& coefList,
449 Int_t code, const RooArgSet* normSet2, const char* rangeName,
450 bool hasWarnedBefore)
451{
452 // Handle trivial passthrough scenario
453 if (code==0) return caller.getVal(normSet2) ;
454
455
456 // WVE needs adaptation for rangeName feature
457 auto* cache = static_cast<CacheElem*>(normIntMgr.getObjByIndex(code-1));
458 if (cache==nullptr) { // revive the (sterilized) cache
459 //cout << "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << ": reviving cache "<< endl;
460 RooArgSet vars;
461 caller.getParameters(nullptr, vars);
462 RooArgSet iset = normIntMgr.selectFromSet2(vars, code-1);
463 RooArgSet nset = normIntMgr.selectFromSet1(vars, code-1);
464 RooArgSet dummy;
465 Int_t code2 = caller.getAnalyticalIntegralWN(iset,dummy,&nset,rangeName);
466 R__ASSERT(code==code2); // must have revived the right (sterilized) slot...
467 cache = static_cast<CacheElem*>(normIntMgr.getObjByIndex(code-1)) ;
468 R__ASSERT(cache!=nullptr);
469 }
470
471 double value(0) ;
472
473 // N funcs, N-1 coefficients
474 double lastCoef(1) ;
475 auto funcIt = funcList.begin();
476 auto funcIntIt = cache->_funcIntList.begin();
477 for (const auto coefArg : coefList) {
479 const auto coef = static_cast<const RooAbsReal*>(coefArg);
480 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
481 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
482
483 double coefVal = coef->getVal(normSet2) ;
484 if (coefVal) {
485 assert(func);
486 if (normSet2 ==nullptr || func->isSelectedComp()) {
488 value += funcInt->getVal()*coefVal ;
489 }
490 lastCoef -= coef->getVal(normSet2) ;
491 }
492 }
493
494 const bool haveLastCoef = funcList.size() == coefList.size();
495
496 if (!haveLastCoef) {
497 // Add last func with correct coefficient
498 const auto func = static_cast<const RooAbsReal*>(*funcIt);
499 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
500 assert(func);
501
502 if (normSet2 ==nullptr || func->isSelectedComp()) {
504 value += funcInt->getVal()*lastCoef ;
505 }
506
507 // Warn about coefficient degeneration
509 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
510 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
511 << 1-lastCoef << std::endl;
512 }
513 }
514
515 double normVal(1) ;
516 if (normSet2 && !normSet2->empty()) {
517 normVal = 0 ;
518
519 // N funcs, N-1 coefficients
520 auto funcNormIter = cache->_funcNormList.begin();
521 for (const auto coefAsArg : coefList) {
522 auto coef = static_cast<RooAbsReal*>(coefAsArg);
523 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
524
525 double coefVal = coef->getVal(normSet2);
526 if (coefVal) {
528 normVal += funcNorm->getVal()*coefVal ;
529 }
530 }
531
532 // Add last func with correct coefficient
533 if (!haveLastCoef) {
534 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
536
537 normVal += funcNorm->getVal()*lastCoef;
538 }
539 }
540
541 return value / normVal;
542}
543
544
545////////////////////////////////////////////////////////////////////////////////
546
548{
549 double n = getNorm(nset) ;
550 if (n<0) {
551 logEvalError("Expected number of events is negative") ;
552 }
553 return n ;
554}
555
556
557////////////////////////////////////////////////////////////////////////////////
558
559std::list<double>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
560{
561 return binBoundaries(_funcList, obs, xlo, xhi);
562}
563
564
565std::list<double> *
566RooRealSumPdf::binBoundaries(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
567{
568 std::unique_ptr<std::list<double>> sumBinB;
569 bool needClean(false);
570
571 // Loop over components pdf
572 for (auto *func : static_range_cast<RooAbsReal *>(funcList)) {
573
574 std::unique_ptr<std::list<double>> funcBinB{func->binBoundaries(obs, xlo, xhi)};
575
576 // Process hint
577 if (funcBinB) {
578 if (!sumBinB) {
579 // If this is the first hint, then just save it
580 sumBinB = std::move(funcBinB);
581 continue;
582 }
583
584 auto newSumBinB = std::make_unique<std::list<double>>(sumBinB->size() + funcBinB->size());
585
586 // Merge hints into temporary array
587 std::merge(funcBinB->begin(), funcBinB->end(), sumBinB->begin(), sumBinB->end(), newSumBinB->begin());
588
589 sumBinB = std::move(newSumBinB);
590 needClean = true;
591 }
592 }
593
594 // Remove consecutive duplicates
595 if (needClean) {
596 sumBinB->erase(std::unique(sumBinB->begin(), sumBinB->end()), sumBinB->end());
597 }
598
599 return sumBinB.release();
600}
601
602
603
604/// Check if all components that depend on `obs` are binned.
606{
607 return isBinnedDistribution(_funcList, obs);
608}
609
610
612{
613 for (auto* func : static_range_cast<RooAbsReal*>(funcList)) {
614
615 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
616 return false ;
617 }
618 }
619
620 return true ;
621}
622
623
624
625////////////////////////////////////////////////////////////////////////////////
626
627std::list<double>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
628{
629 return plotSamplingHint(_funcList, obs, xlo, xhi);
630}
631
632std::list<double> *
633RooRealSumPdf::plotSamplingHint(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
634{
635 std::unique_ptr<std::list<double>> sumHint;
636 bool needClean(false);
637
638 // Loop over components pdf
639 for (const auto elm : funcList) {
640 auto func = static_cast<RooAbsReal *>(elm);
641
642 std::unique_ptr<std::list<double>> funcHint{func->plotSamplingHint(obs, xlo, xhi)};
643
644 // Process hint
645 if (funcHint) {
646 if (!sumHint) {
647
648 // If this is the first hint, then just save it
649 sumHint = std::move(funcHint);
650 continue;
651 }
652
653 auto newSumHint = std::make_unique<std::list<double>>(sumHint->size() + funcHint->size());
654
655 // the lists must be sorted before merging them
656 funcHint->sort();
657 sumHint->sort();
658 // Merge hints into temporary array
659 std::merge(funcHint->begin(), funcHint->end(), sumHint->begin(), sumHint->end(), newSumHint->begin());
660
661 sumHint = std::move(newSumHint);
662 needClean = true;
663 }
664 }
665
666 // Remove consecutive duplicates
667 if (needClean) {
668 sumHint->erase(std::unique(sumHint->begin(), sumHint->end()), sumHint->end());
669 }
670
671 return sumHint.release();
672}
673
674
675
676
677////////////////////////////////////////////////////////////////////////////////
678/// Label OK'ed components of a RooRealSumPdf with cache-and-track
679
684
685
687{
688 for (const auto sarg : funcList) {
689 if (sarg->canNodeBeCached()==Always) {
690 trackNodes.add(*sarg) ;
691 }
692 }
693}
694
695
696////////////////////////////////////////////////////////////////////////////////
697/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
698/// product operator construction
699
700void RooRealSumPdf::printMetaArgs(std::ostream& os) const
701{
703}
704
705
706void RooRealSumPdf::printMetaArgs(RooArgList const& funcList, RooArgList const& coefList, std::ostream& os)
707{
708
709 bool first(true) ;
710
711 if (!coefList.empty()) {
712 auto funcIter = funcList.begin();
713
714 for (const auto coef : coefList) {
715 if (!first) {
716 os << " + " ;
717 } else {
718 first = false ;
719 }
720 const auto func = *(funcIter++);
721 os << coef->GetName() << " * " << func->GetName();
722 }
723
724 if (funcIter != funcList.end()) {
725 os << " + [%] * " << (*funcIter)->GetName() ;
726 }
727 } else {
728
729 for (const auto func : funcList) {
730 if (!first) {
731 os << " + " ;
732 } else {
733 first = false ;
734 }
735 os << func->GetName() ;
736 }
737 }
738
739 os << " " ;
740}
741
742std::unique_ptr<RooAbsArg>
744{
745 if (normSet.empty() || selfNormalized()) {
746 return RooAbsPdf::compileForNormSet({}, ctx);
747 }
748 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
749
750 if (ctx.likelihoodMode() && pdfClone->getAttribute("BinnedLikelihood")) {
751
752 // This has to be done before compiling the servers, such that the
753 // RooBinWidthFunctions know to disable themselves.
754 ctx.setBinnedLikelihoodMode(true);
755
757 ctx.compileServers(*pdfClone, {});
758
759 pdfClone->setAttribute("BinnedLikelihoodActive");
760 // If this is a binned likelihood, we're flagging it in the context.
761 // Then, the RooBinWidthFunctions know that they should not put
762 // themselves in the computation graph. Like this, the pdf values can
763 // directly be interpreted as yields, without multiplying them with the
764 // bin widths again in the NLL. However, the NLL class has to be careful
765 // to only skip the bin with multiplication when there actually were
766 // RooBinWidthFunctions! This is not the case for old workspace before
767 // ROOT 6.26. Therefore, we use the "BinnedLikelihoodActiveYields"
768 // attribute to let the NLL know what it should do.
769 if (ctx.binWidthFuncFlag()) {
770 pdfClone->setAttribute("BinnedLikelihoodActiveYields");
771 }
772 return pdfClone;
773 }
774
775 ctx.compileServers(*pdfClone, {});
776
778 pdfClone->getObservables(&normSet, depList);
779
780 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, depList);
781
782 // The direct servers are this pdf and the normalization integral, which
783 // don't need to be compiled further.
784 for (RooAbsArg *server : newArg->servers()) {
786 }
788 newArg->addOwnedComponents(std::move(pdfClone));
789 return newArg;
790}
791
792std::unique_ptr<RooAbsReal> RooRealSumPdf::createExpectedEventsFunc(const RooArgSet *nset) const
793{
794 if (nset == nullptr)
795 return nullptr;
796 return std::unique_ptr<RooAbsReal>{createIntegral(*nset, *getIntegratorConfig(), normRange())};
797}
#define oocoutW(o, a)
#define coutW(a)
#define oocoutE(o, a)
#define ClassImp(name)
Definition Rtypes.h:382
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#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
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:91
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
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 ...
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.
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
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...
A class to maintain the context for squashing of RooFit models into code.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
void addToCodeBody(RooAbsArg const *klass, std::string const &in)
Adds the input string to the squashed code body.
std::string getTmpVarName() const
Get a unique variable name to be used in the generated code.
std::string buildArg(RooAbsCollection const &x)
Function to save a RooListProxy as an array in the squashed code.
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.
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 std::string translateImpl(RooFit::Detail::CodeSquashContext &ctx, RooAbsArg const *klass, RooArgList const &funcList, RooArgList const &coefList, bool normalize=false)
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
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
const Int_t n
Definition legend1.C:16
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345
static void output()