Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooRealSumPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/** \class RooRealSumPdf
19 \ingroup Roofitcore
20
21
22The class RooRealSumPdf implements a PDF constructed from a sum of functions:
23\f[
24 \mathrm{PDF}(x) = \frac{ \sum_{i=1}^{n-1} \mathrm{coef}_i * \mathrm{func}_i(x) + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \mathrm{func}_n(x) }
25 {\sum_{i=1}^{n-1} \mathrm{coef}_i * \int \mathrm{func}_i(x)dx + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \int \mathrm{func}_n(x) dx }
26\f]
27
28where \f$\mathrm{coef}_i\f$ and \f$\mathrm{func}_i\f$ are RooAbsReal objects, and \f$ x \f$ is the collection of dependents.
29In the present version \f$\mathrm{coef}_i\f$ may not depend on \f$ x \f$, but this limitation could be removed should the need arise.
30
31If the number of coefficients is one less than the number of functions, the PDF is assumed to be normalised. Due to this additional constraint,
32\f$\mathrm{coef}_n\f$ is computed from the other coefficients.
33
34### Extending the PDF
35If an \f$ n^\mathrm{th} \f$ coefficient is provided, the PDF **can** be used as an extended PDF, *i.e.* the total number of events will be measured in addition
36to the fractions of the various functions. **This requires setting the last argument of the constructor to `true`.**
37\note For the RooAddPdf, the extension happens automatically.
38
39### Difference to RooAddPdf / RooRealSumFunc
40- RooAddPdf is a PDF of PDFs, *i.e.* its components need to be normalised and non-negative.
41- RooRealSumPdf is a PDF of functions, *i.e.*, its components can be negative, but their sum cannot be. The normalisation
42 is computed automatically, unless the PDF is extended (see above).
43- RooRealSumFunc is a sum of functions. It is neither normalised, nor need it be positive.
44
45*/
46
47#include "RooRealSumPdf.h"
48
49#include "RooNormalizedPdf.h"
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
62using namespace std;
63
65
67
68////////////////////////////////////////////////////////////////////////////////
69/// Default constructor
70/// coverity[UNINIT_CTOR]
71
72RooRealSumPdf::RooRealSumPdf() : _normIntMgr(this,10)
73{
74
75}
76
77
78
79////////////////////////////////////////////////////////////////////////////////
80/// Constructor with name and title
81
82RooRealSumPdf::RooRealSumPdf(const char *name, const char *title) :
83 RooAbsPdf(name,title),
84 _normIntMgr(this,10),
85 _funcList("!funcList","List of functions",this),
86 _coefList("!coefList","List of coefficients",this),
87 _extended(false),
88 _doFloor(false)
89{
90
91}
92
93
94
95////////////////////////////////////////////////////////////////////////////////
96/// Construct p.d.f consisting of \f$ \mathrm{coef}_1 * \mathrm{func}_1 + (1-\mathrm{coef}_1) * \mathrm{func}_2 \f$.
97/// The input coefficients and functions are allowed to be negative
98/// but the resulting sum is not, which is enforced at runtime.
99
100RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
101 RooAbsReal& func1, RooAbsReal& func2, RooAbsReal& coef1) :
102 RooRealSumPdf(name, title)
103{
104 // Special constructor with two functions and one coefficient
105
106 _funcList.add(func1) ;
107 _funcList.add(func2) ;
108 _coefList.add(coef1) ;
109
110}
111
112
113////////////////////////////////////////////////////////////////////////////////
114/// Constructor for a PDF from a list of functions and coefficients.
115/// It implements
116/// \f[
117/// \sum_i \mathrm{coef}_i \cdot \mathrm{func}_i,
118/// \f]
119/// if \f$ N_\mathrm{coef} = N_\mathrm{func} \f$. With `extended=true`, the coefficients can take any values. With `extended=false`,
120/// there is the danger of getting a degenerate minimisation problem because a PDF has to be normalised, which needs one degree
121/// of freedom less.
122///
123/// A plain (normalised) PDF can therefore be implemented with one less coefficient. RooFit then computes
124/// \f[
125/// \sum_i^{N-1} \mathrm{coef}_i \cdot \mathrm{func}_i + (1 - \sum_i \mathrm{coef}_i ) \cdot \mathrm{func}_N,
126/// \f]
127/// if \f$ N_\mathrm{coef} = N_\mathrm{func} - 1 \f$.
128///
129/// All coefficients and functions are allowed to be negative
130/// but the sum (*i.e.* the PDF) is not, which is enforced at runtime.
131///
132/// \param name Name of the PDF
133/// \param title Title (for plotting)
134/// \param inFuncList List of functions to sum
135/// \param inCoefList List of coefficients
136/// \param extended Interpret as extended PDF (requires equal number of functions and coefficients)
137
138RooRealSumPdf::RooRealSumPdf(const char *name, const char *title,
139 const RooArgList& inFuncList, const RooArgList& inCoefList, bool extended) :
140 RooRealSumPdf(name, title)
141{
142 _extended = extended;
143 RooRealSumPdf::initializeFuncsAndCoefs(*this, inFuncList, inCoefList, _funcList, _coefList);
144}
145
146
148 const RooArgList& inFuncList, const RooArgList& inCoefList,
149 RooArgList& funcList, RooArgList& coefList)
150{
151 const std::string className = caller.ClassName();
152 const std::string constructorName = className + "::" + className;
153
154 if (!(inFuncList.size()==inCoefList.size()+1 || inFuncList.size()==inCoefList.size())) {
155 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName()
156 << ") number of pdfs and coefficients inconsistent, must have Nfunc=Ncoef or Nfunc=Ncoef+1" << std::endl;
157 throw std::invalid_argument(className + ": Number of PDFs and coefficients is inconsistent.");
158 }
159
160 // Constructor with N functions and N or N-1 coefs
161 for (unsigned int i = 0; i < inCoefList.size(); ++i) {
162 const auto& func = inFuncList[i];
163 const auto& coef = inCoefList[i];
164
165 if (!dynamic_cast<const RooAbsReal*>(&coef)) {
166 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") coefficient " << coef.GetName() << " is not of type RooAbsReal, ignored" << std::endl ;
167 continue ;
168 }
169 if (!dynamic_cast<const RooAbsReal*>(&func)) {
170 oocoutW(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") func " << func.GetName() << " is not of type RooAbsReal, ignored" << std::endl ;
171 continue ;
172 }
173 funcList.add(func) ;
174 coefList.add(coef) ;
175 }
176
177 if (inFuncList.size() == inCoefList.size() + 1) {
178 const auto& func = inFuncList[inFuncList.size()-1];
179 if (!dynamic_cast<const RooAbsReal*>(&func)) {
180 oocoutE(&caller, InputArguments) << constructorName << "(" << caller.GetName() << ") last func " << func.GetName() << " is not of type RooAbsReal, fatal error" << std::endl ;
181 throw std::invalid_argument(className + ": Function passed as is not of type RooAbsReal.");
182 }
183 funcList.add(func);
184 }
185
186}
187
188
189
190
191////////////////////////////////////////////////////////////////////////////////
192/// Copy constructor
193
195 RooAbsPdf(other,name),
196 _normIntMgr(other._normIntMgr,this),
197 _funcList("!funcList",this,other._funcList),
198 _coefList("!coefList",this,other._coefList),
199 _extended(other._extended),
200 _doFloor(other._doFloor)
201{
202
203}
204
205
206
207////////////////////////////////////////////////////////////////////////////////
208/// Destructor
209
211{
212
213}
214
215
216
217
218
219////////////////////////////////////////////////////////////////////////////////
220
222{
224}
225
226
228 RooArgList const& funcList,
229 RooArgList const& coefList,
230 bool doFloor,
231 bool & hasWarnedBefore)
232{
233 // Do running sum of coef/func pairs, calculate lastCoef.
234 double value = 0;
235 double sumCoeff = 0.;
236 for (unsigned int i = 0; i < funcList.size(); ++i) {
237 const auto func = static_cast<RooAbsReal*>(&funcList[i]);
238 const auto coef = static_cast<RooAbsReal*>(i < coefList.size() ? &coefList[i] : nullptr);
239 const double coefVal = coef != nullptr ? coef->getVal() : (1. - sumCoeff);
240
241 // Warn about degeneration of last coefficient
242 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
243 if (!hasWarnedBefore) {
244 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
245 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
246 << sumCoeff << ". This means that the PDF is not properly normalised."
247 << " If the PDF was meant to be extended, provide as many coefficients as functions." << std::endl;
248 hasWarnedBefore = true;
249 }
250 // Signal that we are in an undefined region:
251 value = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
252 }
253
254 if (func->isSelectedComp()) {
255 value += func->getVal() * coefVal;
256 }
257
258 sumCoeff += coefVal;
259 }
260
261 // Introduce floor if so requested
262 return value < 0 && doFloor ? 0.0 : value;
263}
264
265
266////////////////////////////////////////////////////////////////////////////////
267/// Calculate the current value
268
270{
272}
273
274
275void RooRealSumPdf::computeBatch(cudaStream_t* /*stream*/, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const {
276 // Do running sum of coef/func pairs, calculate lastCoef.
277 for (unsigned int j = 0; j < nEvents; ++j) {
278 output[j] = 0.0;
279 }
280
281 double sumCoeff = 0.;
282 for (unsigned int i = 0; i < _funcList.size(); ++i) {
283 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
284 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
285 const double coefVal = coef != nullptr ? dataMap.at(coef)[0] : (1. - sumCoeff);
286
287 if (func->isSelectedComp()) {
288 auto funcValues = dataMap.at(func);
289 if(funcValues.size() == 1) {
290 for (unsigned int j = 0; j < nEvents; ++j) {
291 output[j] += funcValues[0] * coefVal;
292 }
293 } else {
294 for (unsigned int j = 0; j < nEvents; ++j) {
295 output[j] += funcValues[j] * coefVal;
296 }
297 }
298 }
299
300 // Warn about degeneration of last coefficient
301 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
302 if (!_haveWarned) {
303 coutW(Eval) << "RooRealSumPdf::evaluateSpan(" << GetName()
304 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
305 << sumCoeff << ". This means that the PDF is not properly normalised. If the PDF was meant to be extended, provide as many coefficients as functions." << endl ;
306 _haveWarned = true;
307 }
308 // Signal that we are in an undefined region by handing back one NaN.
309 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
310 }
311
312 sumCoeff += coefVal;
313 }
314
315 // Introduce floor if so requested
316 if (_doFloor || _doFloorGlobal) {
317 for (unsigned int j = 0; j < nEvents; ++j) {
318 output[j] += std::max(0., output[j]);
319 }
320 }
321}
322
323
325 RooArgList const& funcList, RooArgList const& coefList)
326{
327 bool ret(false) ;
328
329 for (unsigned int i=0; i < coefList.size(); ++i) {
330 const auto& coef = coefList[i];
331 const auto& func = funcList[i];
332
333 if (func.observableOverlaps(nset, coef)) {
334 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
335 << "): ERROR: coefficient " << coef.GetName()
336 << " and FUNC " << func.GetName() << " have one or more observables in common" << std::endl;
337 ret = true ;
338 }
339 if (coef.dependsOn(*nset)) {
340 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
341 << "): ERROR coefficient " << coef.GetName()
342 << " depends on one or more of the following observables" ; nset->Print("1") ;
343 ret = true ;
344 }
345 }
346
347 return ret ;
348}
349
350
351////////////////////////////////////////////////////////////////////////////////
352/// Check if FUNC is valid for given normalization set.
353/// Coefficient and FUNC must be non-overlapping, but func-coefficient
354/// pairs may overlap each other.
355///
356/// In the present implementation, coefficients may not be observables or derive
357/// from observables.
358
360{
361 return checkObservables(*this, nset, _funcList, _coefList);
362}
363
364
365////////////////////////////////////////////////////////////////////////////////
366/// Advertise that all integrals can be handled internally.
367
369 const RooArgSet* normSet2, const char* rangeName) const
370{
371 return getAnalyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, allVars, analVars, normSet2, rangeName);
372}
373
374
376 RooArgList const& funcList, RooArgList const& /*coefList*/,
377 RooArgSet& allVars, RooArgSet& analVars,
378 const RooArgSet* normSet2, const char* rangeName)
379{
380 // Handle trivial no-integration scenario
381 if (allVars.empty()) return 0 ;
382 if (caller.getForceNumInt()) return 0 ;
383
384 // Select subset of allVars that are actual dependents
385 analVars.add(allVars) ;
386 std::unique_ptr<RooArgSet> normSet;
387 if(normSet2) {
388 normSet = std::make_unique<RooArgSet>();
389 caller.getObservables(normSet2, *normSet);
390 }
391
392
393 // Check if this configuration was created before
394 Int_t sterileIdx(-1) ;
395 auto* cache = static_cast<CacheElem*>(normIntMgr.getObj(normSet.get(),&analVars,&sterileIdx,RooNameReg::ptr(rangeName)));
396 if (cache) {
397 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
398 return normIntMgr.lastIndex()+1 ;
399 }
400
401 // Create new cache element
402 cache = new CacheElem ;
403
404 // Make list of function projection and normalization integrals
405 for (const auto elm : funcList) {
406 const auto func = static_cast<RooAbsReal*>(elm);
407
408 std::unique_ptr<RooAbsReal> funcInt{func->createIntegral(analVars,rangeName)};
409 if(auto funcRealInt = dynamic_cast<RooRealIntegral*>(funcInt.get())) funcRealInt->setAllowComponentSelection(true);
410 cache->_funcIntList.addOwned(std::move(funcInt));
411 if (normSet && !normSet->empty()) {
412 cache->_funcNormList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(*normSet)});
413 }
414 }
415
416 // Store cache element
417 Int_t code = normIntMgr.setObj(normSet.get(),&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
418
419 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
420 return code+1 ;
421}
422
423
424
425
426////////////////////////////////////////////////////////////////////////////////
427/// Implement analytical integrations by deferring integration of component
428/// functions to integrators of components.
429
430double RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
431{
432 return analyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, code, normSet2, rangeName, _haveWarned);
433}
434
436 RooArgList const& funcList, RooArgList const& coefList,
437 Int_t code, const RooArgSet* normSet2, const char* rangeName,
438 bool hasWarnedBefore)
439{
440 // Handle trivial passthrough scenario
441 if (code==0) return caller.getVal(normSet2) ;
442
443
444 // WVE needs adaptation for rangeName feature
445 auto* cache = static_cast<CacheElem*>(normIntMgr.getObjByIndex(code-1));
446 if (cache==0) { // revive the (sterilized) cache
447 //cout << "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << ": reviving cache "<< endl;
448 RooArgSet vars;
449 caller.getParameters(nullptr, vars);
450 RooArgSet iset = normIntMgr.selectFromSet2(vars, code-1);
451 RooArgSet nset = normIntMgr.selectFromSet1(vars, code-1);
452 RooArgSet dummy;
453 Int_t code2 = caller.getAnalyticalIntegralWN(iset,dummy,&nset,rangeName);
454 R__ASSERT(code==code2); // must have revived the right (sterilized) slot...
455 cache = (CacheElem*) normIntMgr.getObjByIndex(code-1) ;
456 R__ASSERT(cache!=0);
457 }
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) {
466 assert(funcIt != funcList.end());
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 ==0 || func->isSelectedComp()) {
475 assert(funcInt);
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 ==0 || func->isSelectedComp()) {
491 assert(funcInt);
492 value += funcInt->getVal()*lastCoef ;
493 }
494
495 // Warn about coefficient degeneration
496 if (!hasWarnedBefore && (lastCoef<0 || lastCoef>1)) {
497 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
498 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
499 << 1-lastCoef << endl ;
500 }
501 }
502
503 double normVal(1) ;
504 if (normSet2 && normSet2->getSize()>0) {
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) {
515 assert(funcNorm);
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);
523 assert(funcNorm);
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>* RooRealSumPdf::binBoundaries(RooArgList const& funcList, RooAbsRealLValue& obs, double xlo, double xhi)
554{
555 std::list<double>* sumBinB = nullptr;
556 bool needClean(false) ;
557
558 // Loop over components pdf
559 for (auto * func : static_range_cast<RooAbsReal*>(funcList)) {
560
561 list<double>* funcBinB = func->binBoundaries(obs,xlo,xhi) ;
562
563 // Process hint
564 if (funcBinB) {
565 if (!sumBinB) {
566 // If this is the first hint, then just save it
567 sumBinB = funcBinB ;
568 } else {
569
570 std::list<double>* newSumBinB = new list<double>(sumBinB->size()+funcBinB->size()) ;
571
572 // Merge hints into temporary array
573 merge(funcBinB->begin(),funcBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
574
575 // Copy merged array without duplicates to new sumBinBArrau
576 delete sumBinB ;
577 delete funcBinB ;
578 sumBinB = newSumBinB ;
579 needClean = true ;
580 }
581 }
582 }
583
584 // Remove consecutive duplicates
585 if (needClean) {
586 list<double>::iterator new_end = unique(sumBinB->begin(),sumBinB->end()) ;
587 sumBinB->erase(new_end,sumBinB->end()) ;
588 }
589
590 return sumBinB ;
591}
592
593
594
595/// Check if all components that depend on `obs` are binned.
597{
598 return isBinnedDistribution(_funcList, obs);
599}
600
601
603{
604 for (auto* func : static_range_cast<RooAbsReal*>(funcList)) {
605
606 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
607 return false ;
608 }
609 }
610
611 return true ;
612}
613
614
615
616////////////////////////////////////////////////////////////////////////////////
617
618std::list<double>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
619{
620 return plotSamplingHint(_funcList, obs, xlo, xhi);
621}
622
623std::list<double>* RooRealSumPdf::plotSamplingHint(RooArgList const& funcList, RooAbsRealLValue& obs, double xlo, double xhi)
624{
625 std::list<double>* sumHint = nullptr;
626 bool needClean(false) ;
627
628 // Loop over components pdf
629 for (const auto elm : funcList) {
630 auto func = static_cast<RooAbsReal*>(elm);
631
632 list<double>* funcHint = func->plotSamplingHint(obs,xlo,xhi) ;
633
634 // Process hint
635 if (funcHint) {
636 if (!sumHint) {
637
638 // If this is the first hint, then just save it
639 sumHint = funcHint ;
640
641 } else {
642
643 auto* newSumHint = new std::list<double>(sumHint->size()+funcHint->size()) ;
644
645 // the lists must be sorted before merging them
646 funcHint->sort();
647 sumHint->sort();
648 // Merge hints into temporary array
649 merge(funcHint->begin(),funcHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
650
651 // Copy merged array without duplicates to new sumHintArrau
652 delete sumHint ;
653 sumHint = newSumHint ;
654 needClean = true ;
655 }
656 }
657 }
658
659 // Remove consecutive duplicates
660 if (needClean) {
661 sumHint->erase(std::unique(sumHint->begin(),sumHint->end()), sumHint->end()) ;
662 }
663
664 return sumHint ;
665}
666
667
668
669
670////////////////////////////////////////////////////////////////////////////////
671/// Label OK'ed components of a RooRealSumPdf with cache-and-track
672
674{
675 setCacheAndTrackHints(_funcList, trackNodes);
676}
677
678
680{
681 for (const auto sarg : funcList) {
682 if (sarg->canNodeBeCached()==Always) {
683 trackNodes.add(*sarg) ;
684 //cout << "tracking node RealSumPdf component " << sarg->ClassName() << "::" << sarg->GetName() << endl ;
685 }
686 }
687}
688
689
690////////////////////////////////////////////////////////////////////////////////
691/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
692/// product operator construction
693
694void RooRealSumPdf::printMetaArgs(ostream& os) const
695{
697}
698
699
700void RooRealSumPdf::printMetaArgs(RooArgList const& funcList, RooArgList const& coefList, ostream& os)
701{
702
703 bool first(true) ;
704
705 if (!coefList.empty()) {
706 auto funcIter = funcList.begin();
707
708 for (const auto coef : coefList) {
709 if (!first) {
710 os << " + " ;
711 } else {
712 first = false ;
713 }
714 const auto func = *(funcIter++);
715 os << coef->GetName() << " * " << func->GetName();
716 }
717
718 if (funcIter != funcList.end()) {
719 os << " + [%] * " << (*funcIter)->GetName() ;
720 }
721 } else {
722
723 for (const auto func : funcList) {
724 if (!first) {
725 os << " + " ;
726 } else {
727 first = false ;
728 }
729 os << func->GetName() ;
730 }
731 }
732
733 os << " " ;
734}
735
736std::unique_ptr<RooAbsArg>
738{
739 if (normSet.empty() || selfNormalized()) {
740 return RooAbsPdf::compileForNormSet({}, ctx);
741 }
742 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
743
744 if (ctx.likelihoodMode() && pdfClone->getAttribute("BinnedLikelihood")) {
745
746 // This has to be done before compiling the servers, such that the
747 // RooBinWidthFunctions know to disable themselves.
748 ctx.setBinnedLikelihoodMode(true);
749
750 ctx.markAsCompiled(*pdfClone);
751 ctx.compileServers(*pdfClone, {});
752
753 pdfClone->setAttribute("BinnedLikelihoodActive");
754 // If this is a binned likelihood, we're flagging it in the context.
755 // Then, the RooBinWidthFunctions know that they should not put
756 // themselves in the computation graph. Like this, the pdf values can
757 // directly be interpreted as yields, without multiplying them with the
758 // bin widths again in the NLL. However, the NLL class has to be careful
759 // to only skip the bin with multiplication when there actually were
760 // RooBinWidthFunctions! This is not the case for old workspace before
761 // ROOT 6.26. Therefore, we use the "BinnedLikelihoodActiveYields"
762 // attribute to let the NLL know what it should do.
763 if (ctx.binWidthFuncFlag()) {
764 pdfClone->setAttribute("BinnedLikelihoodActiveYields");
765 }
766 return pdfClone;
767 }
768
769 ctx.compileServers(*pdfClone, {});
770
771 RooArgSet depList;
772 pdfClone->getObservables(&normSet, depList);
773
774 auto newArg = std::make_unique<RooNormalizedPdf>(*pdfClone, depList);
775
776 // The direct servers are this pdf and the normalization integral, which
777 // don't need to be compiled further.
778 for (RooAbsArg *server : newArg->servers()) {
779 ctx.markAsCompiled(*server);
780 }
781 ctx.markAsCompiled(*newArg);
782 newArg->addOwnedComponents(std::move(pdfClone));
783 return newArg;
784}
785
786std::unique_ptr<RooAbsReal> RooRealSumPdf::createExpectedEventsFunc(const RooArgSet *nset) const
787{
788 if (nset == nullptr)
789 return nullptr;
790 return std::unique_ptr<RooAbsReal>{createIntegral(*nset, *getIntegratorConfig(), normRange())};
791}
#define oocoutW(o, a)
#define coutW(a)
#define oocoutE(o, a)
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:117
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
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
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:86
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
Int_t getSize() const
Return the number of elements in the collection.
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.
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:255
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
@ CanBeExtended
Definition RooAbsPdf.h:272
@ CanNotBeExtended
Definition RooAbsPdf.h:272
const char * normRange() const
Definition RooAbsPdf.h:310
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
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 ...
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
bool getForceNumInt() const
Definition RooAbsReal.h:178
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...
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:55
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 contatining the objects that are both in the cached set 1.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 2.
Int_t lastIndex() const
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)
RooSpan< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:21
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
The class RooRealSumPdf 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.
const RooArgList & funcList() const
~RooRealSumPdf() override
Destructor.
RooObjCacheManager _normIntMgr
! The integration cache manager
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const override
Base function for computing multiple values of a RooAbsReal.
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...
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 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:207
const Int_t n
Definition legend1.C:16
Definition first.py:1
__roodevice__ static __roohost__ double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output()