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
62
64
65////////////////////////////////////////////////////////////////////////////////
66/// Default constructor
67/// coverity[UNINIT_CTOR]
68
69RooRealSumPdf::RooRealSumPdf() : _normIntMgr(this, 10) {}
70
71////////////////////////////////////////////////////////////////////////////////
72/// Constructor with name and title
73
74RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
75 RooAbsPdf(name,title),
76 _normIntMgr(this,10),
77 _funcList("!funcList","List of functions",this),
78 _coefList("!coefList","List of coefficients",this),
79 _extended(false)
80{
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Construct p.d.f consisting of \f$ \mathrm{coef}_1 * \mathrm{func}_1 + (1-\mathrm{coef}_1) * \mathrm{func}_2 \f$.
85/// The input coefficients and functions are allowed to be negative
86/// but the resulting sum is not, which is enforced at runtime.
87
88RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
90 RooRealSumPdf(name, title)
91{
92 // Special constructor with two functions and one coefficient
93
97
98}
99
100
101////////////////////////////////////////////////////////////////////////////////
102/// Constructor for a PDF from a list of functions and coefficients.
103/// It implements
104/// \f[
105/// \sum_i \mathrm{coef}_i \cdot \mathrm{func}_i,
106/// \f]
107/// if \f$ N_\mathrm{coef} = N_\mathrm{func} \f$. With `extended=true`, the coefficients can take any values. With `extended=false`,
108/// there is the danger of getting a degenerate minimisation problem because a PDF has to be normalised, which needs one degree
109/// of freedom less.
110///
111/// A plain (normalised) PDF can therefore be implemented with one less coefficient. RooFit then computes
112/// \f[
113/// \sum_i^{N-1} \mathrm{coef}_i \cdot \mathrm{func}_i + (1 - \sum_i \mathrm{coef}_i ) \cdot \mathrm{func}_N,
114/// \f]
115/// if \f$ N_\mathrm{coef} = N_\mathrm{func} - 1 \f$.
116///
117/// All coefficients and functions are allowed to be negative
118/// but the sum (*i.e.* the PDF) is not, which is enforced at runtime.
119///
120/// \param name Name of the PDF
121/// \param title Title (for plotting)
122/// \param inFuncList List of functions to sum
123/// \param inCoefList List of coefficients
124/// \param extended Interpret as extended PDF (requires equal number of functions and coefficients)
125
134
135
138 RooArgList& funcList, RooArgList& coefList)
139{
140 const std::string className = caller.ClassName();
141 const std::string constructorName = className + "::" + className;
142
143 if (!(inFuncList.size()==inCoefList.size()+1 || inFuncList.size()==inCoefList.size())) {
144 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName()
145 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << std::endl;
146 throw std::invalid_argument(className + ": Number of PDFs and coefficients is inconsistent.");
147 }
148
149 // Constructor with N functions and N or N-1 coefs
150 for (unsigned int i = 0; i < inCoefList.size(); ++i) {
151 const auto& func = inFuncList[i];
152 const auto& coef = inCoefList[i];
153
154 if (!dynamic_cast<const RooAbsReal*>(&coef)) {
155 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") coefficient " << coef.GetName() << " is not of type RooAbsReal, ignored" << std::endl ;
156 continue ;
157 }
158 if (!dynamic_cast<const RooAbsReal*>(&func)) {
159 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") func " << func.GetName() << " is not of type RooAbsReal, ignored" << std::endl ;
160 continue ;
161 }
162 funcList.add(func) ;
163 coefList.add(coef) ;
164 }
165
166 if (inFuncList.size() == inCoefList.size() + 1) {
167 const auto& func = inFuncList[inFuncList.size()-1];
168 if (!dynamic_cast<const RooAbsReal*>(&func)) {
169 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") last func " << func.GetName() << " is not of type RooAbsReal, fatal error" << std::endl ;
170 throw std::invalid_argument(className + ": Function passed as is not of type RooAbsReal.");
171 }
172 funcList.add(func);
173 }
174
175}
176
177
178
179
180////////////////////////////////////////////////////////////////////////////////
181/// Copy constructor
182
185 _normIntMgr(other._normIntMgr,this),
186 _funcList("!funcList",this,other._funcList),
187 _coefList("!coefList",this,other._coefList),
188 _extended(other._extended),
189 _doFloor(other._doFloor)
190{
191
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// Collect the list of functions to be integrated from the cache
196
198{
199 return getCacheElem(*this, _normIntMgr, code, rangeName)->_funcIntList;
200}
201
203 Int_t code, const char *rangeName)
204{
205 return getCacheElem(caller, normIntMgr, code, rangeName)->_funcIntList;
206}
207
208////////////////////////////////////////////////////////////////////////////////
209
214
215////////////////////////////////////////////////////////////////////////////////
216
218 RooArgList const& funcList,
219 RooArgList const& coefList,
220 bool doFloor,
221 bool & hasWarnedBefore)
222{
223 // Do running sum of coef/func pairs, calculate lastCoef.
224 double value = 0;
225 double sumCoeff = 0.;
226 for (unsigned int i = 0; i < funcList.size(); ++i) {
227 const auto func = static_cast<RooAbsReal*>(&funcList[i]);
228 const auto coef = static_cast<RooAbsReal*>(i < coefList.size() ? &coefList[i] : nullptr);
229 const double coefVal = coef != nullptr ? coef->getVal() : (1. - sumCoeff);
230
231 // Warn about degeneration of last coefficient
232 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
233 if (!hasWarnedBefore) {
234 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
235 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
236 << sumCoeff << ". This means that the PDF is not properly normalised."
237 << " If the PDF was meant to be extended, provide as many coefficients as functions." << std::endl;
238 hasWarnedBefore = 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 return value < 0 && doFloor ? 0.0 : value;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// Calculate the current value
257
259{
261}
262
263
265{
266 std::span<double> output = ctx.output();
267 std::size_t nEvents = output.size();
268
269 // Do running sum of coef/func pairs, calculate lastCoef.
270 std::fill(output.begin(), output.end(), 0.);
271
272 double sumCoeff = 0.;
273 for (unsigned int i = 0; i < _funcList.size(); ++i) {
274 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
275 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
276 const double coefVal = coef != nullptr ? ctx.at(coef)[0] : (1. - sumCoeff);
277
278 if (func->isSelectedComp()) {
279 auto funcValues = ctx.at(func);
280 if(funcValues.size() == 1) {
281 for (unsigned int j = 0; j < nEvents; ++j) {
282 output[j] += funcValues[0] * coefVal;
283 }
284 } else {
285 for (unsigned int j = 0; j < nEvents; ++j) {
286 output[j] += funcValues[j] * coefVal;
287 }
288 }
289 }
290
291 // Warn about degeneration of last coefficient
292 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
293 if (!_haveWarned) {
294 coutW(Eval) << "RooRealSumPdf::doEval(" << GetName()
295 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
296 << 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;
297 _haveWarned = true;
298 }
299 // Signal that we are in an undefined region by handing back one NaN.
300 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
301 }
302
303 sumCoeff += coefVal;
304 }
305
306 // Introduce floor if so requested
307 if (_doFloor || _doFloorGlobal) {
308 for (unsigned int j = 0; j < nEvents; ++j) {
309 output[j] = std::max(0., output[j]);
310 }
311 }
312}
313
315 RooArgList const& funcList, RooArgList const& coefList)
316{
317 bool ret(false) ;
318
319 for (unsigned int i=0; i < coefList.size(); ++i) {
320 const auto& coef = coefList[i];
321 const auto& func = funcList[i];
322
323 if (func.observableOverlaps(nset, coef)) {
324 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
325 << "): ERROR: coefficient " << coef.GetName()
326 << " and FUNC " << func.GetName() << " have one or more observables in common" << std::endl;
327 ret = true ;
328 }
329 if (coef.dependsOn(*nset)) {
330 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
331 << "): ERROR coefficient " << coef.GetName()
332 << " depends on one or more of the following observables" ; nset->Print("1") ;
333 ret = true ;
334 }
335 }
336
337 return ret ;
338}
339
340
341////////////////////////////////////////////////////////////////////////////////
342/// Check if FUNC is valid for given normalization set.
343/// Coefficient and FUNC must be non-overlapping, but func-coefficient
344/// pairs may overlap each other.
345///
346/// In the present implementation, coefficients may not be observables or derive
347/// from observables.
348
350{
351 return checkObservables(*this, nset, _funcList, _coefList);
352}
353
354
355////////////////////////////////////////////////////////////////////////////////
356/// Advertise that all integrals can be handled internally.
357
363
364
366 RooArgList const& funcList, RooArgList const& /*coefList*/,
367 RooArgSet& allVars, RooArgSet& analVars,
368 const RooArgSet* normSet2, const char* rangeName)
369{
370 // Handle trivial no-integration scenario
371 if (allVars.empty()) return 0 ;
372 if (caller.getForceNumInt()) return 0 ;
373
374 // Select subset of allVars that are actual dependents
375 analVars.add(allVars) ;
376 std::unique_ptr<RooArgSet> normSet;
377 if(normSet2) {
378 normSet = std::make_unique<RooArgSet>();
379 caller.getObservables(normSet2, *normSet);
380 }
381
382
383 // Check if this configuration was created before
384 Int_t sterileIdx(-1) ;
385 auto* cache = static_cast<CacheElem*>(normIntMgr.getObj(normSet.get(),&analVars,&sterileIdx,RooNameReg::ptr(rangeName)));
386 if (cache) {
387 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << std::endl;
388 return normIntMgr.lastIndex()+1 ;
389 }
390
391 // Create new cache element
392 cache = new CacheElem ;
393
394 // Make list of function projection and normalization integrals
395 for (const auto elm : funcList) {
396 const auto func = static_cast<RooAbsReal*>(elm);
397
398 std::unique_ptr<RooAbsReal> funcInt{func->createIntegral(analVars,rangeName)};
399 if(auto funcRealInt = dynamic_cast<RooRealIntegral*>(funcInt.get())) funcRealInt->setAllowComponentSelection(true);
400 cache->_funcIntList.addOwned(std::move(funcInt));
401 if (normSet && !normSet->empty()) {
402 cache->_funcNormList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(*normSet)});
403 }
404 }
405
406 // Store cache element
408
409 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << std::endl;
410 return code+1 ;
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// Collect cache elements for given code to use in codegen
415
418{
419 // WVE needs adaptation for rangeName feature
420 auto *cache = static_cast<CacheElem *>(normIntMgr.getObjByIndex(code - 1));
421 if (cache == nullptr) { // revive the (sterilized) cache
422 // std::cout <<
423 // "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>")
424 // << ": reviving cache "<< std::endl;
425 RooArgSet vars;
426 caller.getParameters(nullptr, vars);
427 RooArgSet iset = normIntMgr.selectFromSet2(vars, code - 1);
428 RooArgSet nset = normIntMgr.selectFromSet1(vars, code - 1);
429 RooArgSet dummy;
430 Int_t code2 = caller.getAnalyticalIntegralWN(iset, dummy, &nset, rangeName);
431 R__ASSERT(code == code2); // must have revived the right (sterilized) slot...
432 cache = static_cast<CacheElem *>(normIntMgr.getObjByIndex(code - 1));
433 R__ASSERT(cache != nullptr);
434 }
435
436 return cache;
437}
438
439////////////////////////////////////////////////////////////////////////////////
440/// Implement analytical integrations by deferring integration of component
441/// functions to integrators of components.
442
447
449 RooArgList const& funcList, RooArgList const& coefList,
450 Int_t code, const RooArgSet* normSet2, const char* rangeName,
451 bool hasWarnedBefore)
452{
453 // Handle trivial passthrough scenario
454 if (code==0) return caller.getVal(normSet2) ;
455
456 // Retrieve cache element
457 const CacheElem *cache = getCacheElem(caller, normIntMgr, code, rangeName);
458
459 double value(0) ;
460
461 // N funcs, N-1 coefficients
462 double lastCoef(1) ;
463 auto funcIt = funcList.begin();
464 auto funcIntIt = cache->_funcIntList.begin();
465 for (const auto coefArg : coefList) {
467 const auto coef = static_cast<const RooAbsReal*>(coefArg);
468 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
469 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
470
471 double coefVal = coef->getVal(normSet2) ;
472 if (coefVal) {
473 assert(func);
474 if (normSet2 ==nullptr || func->isSelectedComp()) {
476 value += funcInt->getVal()*coefVal ;
477 }
478 lastCoef -= coef->getVal(normSet2) ;
479 }
480 }
481
482 const bool haveLastCoef = funcList.size() == coefList.size();
483
484 if (!haveLastCoef) {
485 // Add last func with correct coefficient
486 const auto func = static_cast<const RooAbsReal*>(*funcIt);
487 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
488 assert(func);
489
490 if (normSet2 ==nullptr || func->isSelectedComp()) {
492 value += funcInt->getVal()*lastCoef ;
493 }
494
495 // Warn about coefficient degeneration
497 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
498 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
499 << 1-lastCoef << std::endl;
500 }
501 }
502
503 double normVal(1) ;
504 if (normSet2 && !normSet2->empty()) {
505 normVal = 0 ;
506
507 // N funcs, N-1 coefficients
508 auto funcNormIter = cache->_funcNormList.begin();
509 for (const auto coefAsArg : coefList) {
510 auto coef = static_cast<RooAbsReal*>(coefAsArg);
511 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
512
513 double coefVal = coef->getVal(normSet2);
514 if (coefVal) {
516 normVal += funcNorm->getVal()*coefVal ;
517 }
518 }
519
520 // Add last func with correct coefficient
521 if (!haveLastCoef) {
522 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
524
525 normVal += funcNorm->getVal()*lastCoef;
526 }
527 }
528
529 return value / normVal;
530}
531
532
533////////////////////////////////////////////////////////////////////////////////
534
536{
537 double n = getNorm(nset) ;
538 if (n<0) {
539 logEvalError("Expected number of events is negative") ;
540 }
541 return n ;
542}
543
544
545////////////////////////////////////////////////////////////////////////////////
546
547std::list<double>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
548{
549 return binBoundaries(_funcList, obs, xlo, xhi);
550}
551
552
553std::list<double> *
554RooRealSumPdf::binBoundaries(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
555{
556 std::unique_ptr<std::list<double>> sumBinB;
557 bool needClean(false);
558
559 // Loop over components pdf
560 for (auto *func : static_range_cast<RooAbsReal *>(funcList)) {
561
562 std::unique_ptr<std::list<double>> funcBinB{func->binBoundaries(obs, xlo, xhi)};
563
564 // Process hint
565 if (funcBinB) {
566 if (!sumBinB) {
567 // If this is the first hint, then just save it
568 sumBinB = std::move(funcBinB);
569 continue;
570 }
571
572 auto newSumBinB = std::make_unique<std::list<double>>(sumBinB->size() + funcBinB->size());
573
574 // Merge hints into temporary array
575 std::merge(funcBinB->begin(), funcBinB->end(), sumBinB->begin(), sumBinB->end(), newSumBinB->begin());
576
577 sumBinB = std::move(newSumBinB);
578 needClean = true;
579 }
580 }
581
582 // Remove consecutive duplicates
583 if (needClean) {
584 sumBinB->erase(std::unique(sumBinB->begin(), sumBinB->end()), sumBinB->end());
585 }
586
587 return sumBinB.release();
588}
589
590
591
592/// Check if all components that depend on `obs` are binned.
594{
595 return isBinnedDistribution(_funcList, obs);
596}
597
598
600{
601 for (auto* func : static_range_cast<RooAbsReal*>(funcList)) {
602
603 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
604 return false ;
605 }
606 }
607
608 return true ;
609}
610
611
612
613////////////////////////////////////////////////////////////////////////////////
614
615std::list<double>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
616{
617 return plotSamplingHint(_funcList, obs, xlo, xhi);
618}
619
620std::list<double> *
621RooRealSumPdf::plotSamplingHint(RooArgList const &funcList, RooAbsRealLValue &obs, double xlo, double xhi)
622{
623 std::unique_ptr<std::list<double>> sumHint;
624 bool needClean(false);
625
626 // Loop over components pdf
627 for (const auto elm : funcList) {
628 auto func = static_cast<RooAbsReal *>(elm);
629
630 std::unique_ptr<std::list<double>> funcHint{func->plotSamplingHint(obs, xlo, xhi)};
631
632 // Process hint
633 if (funcHint) {
634 if (!sumHint) {
635
636 // If this is the first hint, then just save it
637 sumHint = std::move(funcHint);
638 continue;
639 }
640
641 auto newSumHint = std::make_unique<std::list<double>>(sumHint->size() + funcHint->size());
642
643 // the lists must be sorted before merging them
644 funcHint->sort();
645 sumHint->sort();
646 // Merge hints into temporary array
647 std::merge(funcHint->begin(), funcHint->end(), sumHint->begin(), sumHint->end(), newSumHint->begin());
648
649 sumHint = std::move(newSumHint);
650 needClean = true;
651 }
652 }
653
654 // Remove consecutive duplicates
655 if (needClean) {
656 sumHint->erase(std::unique(sumHint->begin(), sumHint->end()), sumHint->end());
657 }
658
659 return sumHint.release();
660}
661
662
663
664
665////////////////////////////////////////////////////////////////////////////////
666/// Label OK'ed components of a RooRealSumPdf with cache-and-track
667
672
673
675{
676 for (const auto sarg : funcList) {
677 if (sarg->canNodeBeCached()==Always) {
678 trackNodes.add(*sarg) ;
679 }
680 }
681}
682
683
684////////////////////////////////////////////////////////////////////////////////
685/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
686/// product operator construction
687
688void RooRealSumPdf::printMetaArgs(std::ostream& os) const
689{
691}
692
693
694void RooRealSumPdf::printMetaArgs(RooArgList const& funcList, RooArgList const& coefList, std::ostream& os)
695{
696
697 bool first(true) ;
698
699 if (!coefList.empty()) {
700 auto funcIter = funcList.begin();
701
702 for (const auto coef : coefList) {
703 if (!first) {
704 os << " + " ;
705 } else {
706 first = false ;
707 }
708 const auto func = *(funcIter++);
709 os << coef->GetName() << " * " << func->GetName();
710 }
711
712 if (funcIter != funcList.end()) {
713 os << " + [%] * " << (*funcIter)->GetName() ;
714 }
715 } else {
716
717 for (const auto func : funcList) {
718 if (!first) {
719 os << " + " ;
720 } else {
721 first = false ;
722 }
723 os << func->GetName() ;
724 }
725 }
726
727 os << " " ;
728}
729
730std::unique_ptr<RooAbsArg>
732{
733 if (normSet.empty() || selfNormalized()) {
734 return RooAbsPdf::compileForNormSet({}, ctx);
735 }
736 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
737
738 if (ctx.likelihoodMode() && pdfClone->getAttribute("BinnedLikelihood")) {
739
740 // This has to be done before compiling the servers, such that the
741 // RooBinWidthFunctions know to disable themselves.
742 ctx.setBinnedLikelihoodMode(true);
743
745 ctx.compileServers(*pdfClone, {});
746
747 pdfClone->setAttribute("BinnedLikelihoodActive");
748 // If this is a binned likelihood, we're flagging it in the context.
749 // Then, the RooBinWidthFunctions know that they should not put
750 // themselves in the computation graph. Like this, the pdf values can
751 // directly be interpreted as yields, without multiplying them with the
752 // bin widths again in the NLL. However, the NLL class has to be careful
753 // to only skip the bin with multiplication when there actually were
754 // RooBinWidthFunctions! This is not the case for old workspace before
755 // ROOT 6.26. Therefore, we use the "BinnedLikelihoodActiveYields"
756 // attribute to let the NLL know what it should do.
757 if (ctx.binWidthFuncFlag()) {
758 pdfClone->setAttribute("BinnedLikelihoodActiveYields");
759 }
760 return pdfClone;
761 }
762
763 ctx.compileServers(*pdfClone, {});
764
766 pdfClone->getObservables(&normSet, depList);
767
768 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, depList);
769
770 // The direct servers are this pdf and the normalization integral, which
771 // don't need to be compiled further.
772 for (RooAbsArg *server : newArg->servers()) {
774 }
776 newArg->addOwnedComponents(std::move(pdfClone));
777 return newArg;
778}
779
780std::unique_ptr<RooAbsReal> RooRealSumPdf::createExpectedEventsFunc(const RooArgSet *nset) const
781{
782 if (nset == nullptr)
783 return nullptr;
784 return std::unique_ptr<RooAbsReal>{createIntegral(*nset, *getIntegratorConfig(), normRange())};
785}
#define oocoutW(o, a)
#define coutW(a)
#define oocoutE(o, a)
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:77
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.
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
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.
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
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:49
const Int_t n
Definition legend1.C:16
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output()