Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooRealSumPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/** \class RooRealSumPdf
19 \ingroup Roofitcore
20
21
22Implements a PDF constructed from a sum of functions:
23\f[
24 \mathrm{PDF}(x) = \frac{ \sum_{i=1}^{n-1} \mathrm{coef}_i * \mathrm{func}_i(x) + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \mathrm{func}_n(x) }
25 {\sum_{i=1}^{n-1} \mathrm{coef}_i * \int \mathrm{func}_i(x)dx + \left[ 1 - \sum_{i=1}^{n-1} \mathrm{coef}_i \right] * \int \mathrm{func}_n(x) dx }
26\f]
27
28where \f$\mathrm{coef}_i\f$ and \f$\mathrm{func}_i\f$ are RooAbsReal objects, and \f$ x \f$ is the collection of dependents.
29In the present version \f$\mathrm{coef}_i\f$ may not depend on \f$ x \f$, but this limitation could be removed should the need arise.
30
31If the number of coefficients is one less than the number of functions, the PDF is assumed to be normalised. Due to this additional constraint,
32\f$\mathrm{coef}_n\f$ is computed from the other coefficients.
33
34### Extending the PDF
35If an \f$ n^\mathrm{th} \f$ coefficient is provided, the PDF **can** be used as an extended PDF, *i.e.* the total number of events will be measured in addition
36to the fractions of the various functions. **This requires setting the last argument of the constructor to `true`.**
37\note For the RooAddPdf, the extension happens automatically.
38
39### Difference to RooAddPdf / RooRealSumFunc
40- RooAddPdf is a PDF of PDFs, *i.e.* its components need to be normalised and non-negative.
41- RooRealSumPdf is a PDF of functions, *i.e.*, its components can be negative, but their sum cannot be. The normalisation
42 is computed automatically, unless the PDF is extended (see above).
43- RooRealSumFunc is a sum of functions. It is neither normalised, nor need it be positive.
44
45*/
46
47#include "RooRealSumPdf.h"
48
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/// Calculate the current value
267
269{
271}
272
273
274void RooRealSumPdf::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const {
275 // Do running sum of coef/func pairs, calculate lastCoef.
276 for (unsigned int j = 0; j < nEvents; ++j) {
277 output[j] = 0.0;
278 }
279
280 double sumCoeff = 0.;
281 for (unsigned int i = 0; i < _funcList.size(); ++i) {
282 const auto func = static_cast<RooAbsReal*>(&_funcList[i]);
283 const auto coef = static_cast<RooAbsReal*>(i < _coefList.size() ? &_coefList[i] : nullptr);
284 const double coefVal = coef != nullptr ? dataMap.at(coef)[0] : (1. - sumCoeff);
285
286 if (func->isSelectedComp()) {
287 auto funcValues = dataMap.at(func);
288 if(funcValues.size() == 1) {
289 for (unsigned int j = 0; j < nEvents; ++j) {
290 output[j] += funcValues[0] * coefVal;
291 }
292 } else {
293 for (unsigned int j = 0; j < nEvents; ++j) {
294 output[j] += funcValues[j] * coefVal;
295 }
296 }
297 }
298
299 // Warn about degeneration of last coefficient
300 if (coef == nullptr && (coefVal < 0 || coefVal > 1.)) {
301 if (!_haveWarned) {
302 coutW(Eval) << "RooRealSumPdf::computeBatch(" << GetName()
303 << ") WARNING: sum of FUNC coefficients not in range [0-1], value="
304 << 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 ;
305 _haveWarned = true;
306 }
307 // Signal that we are in an undefined region by handing back one NaN.
308 output[0] = RooNaNPacker::packFloatIntoNaN(100.f * (coefVal < 0. ? -coefVal : coefVal - 1.));
309 }
310
311 sumCoeff += coefVal;
312 }
313
314 // Introduce floor if so requested
315 if (_doFloor || _doFloorGlobal) {
316 for (unsigned int j = 0; j < nEvents; ++j) {
317 output[j] += std::max(0., output[j]);
318 }
319 }
320}
321
323 RooArgList const &funcList, RooArgList const &coefList)
324{
325 bool noLastCoeff = funcList.size() != coefList.size();
326 std::string sum;
327 std::string lastCoeff;
328
329 if (funcList.size() < 3) {
330 lastCoeff = "(1";
331
332 std::size_t i;
333 for (i = 0; i < coefList.size(); ++i) {
334 auto coeff = ctx.getResult(coefList[i]);
335 sum += "(" + ctx.getResult(funcList[i]) + " * " + coeff + ") +";
336 if (noLastCoeff)
337 lastCoeff += " - " + coeff;
338 }
339 lastCoeff += ")";
340
341 if (noLastCoeff) {
342 sum += "(" + ctx.getResult(funcList[i]) + " * " + lastCoeff + ")";
343 } else {
344 sum.pop_back();
345 }
346
347 } else {
348 std::string const &funcName = ctx.buildArg(funcList);
349 std::string const &coeffName = ctx.buildArg(coefList);
350 std::string const &coeffSize = std::to_string(coefList.size());
351 std::string const &className = klass->GetName();
352
353 sum = ctx.getTmpVarName();
354 lastCoeff = ctx.getTmpVarName();
355 ctx.addToCodeBody(klass, "double " + sum + " = 0, " + lastCoeff + "= 0;\n");
356
357 std::string iterator = "i_" + className;
358 std::string subscriptExpr = "[" + iterator + "]";
359
360 std::string code = "for(int " + iterator + " = 0; " + iterator + " < " + coeffSize + "; " + iterator + "++) {\n" +
361 sum + " += " + funcName + subscriptExpr + " * " + coeffName + subscriptExpr + ";\n";
362 if (noLastCoeff)
363 code += lastCoeff + " += " + coeffName + subscriptExpr + ";\n";
364 code += "}\n";
365
366 if (noLastCoeff)
367 code += sum + " += " + funcName + "[" + coeffSize + "]" + " * (1 - " + lastCoeff + ");\n";
368 ctx.addToCodeBody(klass, code);
369 }
370
371 ctx.addResult(klass, sum);
372}
373
375{
376 translateImpl(ctx, this, _funcList, _coefList);
377}
378
380 RooArgList const& funcList, RooArgList const& coefList)
381{
382 bool ret(false) ;
383
384 for (unsigned int i=0; i < coefList.size(); ++i) {
385 const auto& coef = coefList[i];
386 const auto& func = funcList[i];
387
388 if (func.observableOverlaps(nset, coef)) {
389 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
390 << "): ERROR: coefficient " << coef.GetName()
391 << " and FUNC " << func.GetName() << " have one or more observables in common" << std::endl;
392 ret = true ;
393 }
394 if (coef.dependsOn(*nset)) {
395 oocoutE(&caller, InputArguments) << caller.ClassName() << "::checkObservables(" << caller.GetName()
396 << "): ERROR coefficient " << coef.GetName()
397 << " depends on one or more of the following observables" ; nset->Print("1") ;
398 ret = true ;
399 }
400 }
401
402 return ret ;
403}
404
405
406////////////////////////////////////////////////////////////////////////////////
407/// Check if FUNC is valid for given normalization set.
408/// Coefficient and FUNC must be non-overlapping, but func-coefficient
409/// pairs may overlap each other.
410///
411/// In the present implementation, coefficients may not be observables or derive
412/// from observables.
413
415{
416 return checkObservables(*this, nset, _funcList, _coefList);
417}
418
419
420////////////////////////////////////////////////////////////////////////////////
421/// Advertise that all integrals can be handled internally.
422
424 const RooArgSet* normSet2, const char* rangeName) const
425{
426 return getAnalyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, allVars, analVars, normSet2, rangeName);
427}
428
429
431 RooArgList const& funcList, RooArgList const& /*coefList*/,
432 RooArgSet& allVars, RooArgSet& analVars,
433 const RooArgSet* normSet2, const char* rangeName)
434{
435 // Handle trivial no-integration scenario
436 if (allVars.empty()) return 0 ;
437 if (caller.getForceNumInt()) return 0 ;
438
439 // Select subset of allVars that are actual dependents
440 analVars.add(allVars) ;
441 std::unique_ptr<RooArgSet> normSet;
442 if(normSet2) {
443 normSet = std::make_unique<RooArgSet>();
444 caller.getObservables(normSet2, *normSet);
445 }
446
447
448 // Check if this configuration was created before
449 Int_t sterileIdx(-1) ;
450 auto* cache = static_cast<CacheElem*>(normIntMgr.getObj(normSet.get(),&analVars,&sterileIdx,RooNameReg::ptr(rangeName)));
451 if (cache) {
452 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << _normIntMgr.lastIndex()+1 << " (cached)" << endl;
453 return normIntMgr.lastIndex()+1 ;
454 }
455
456 // Create new cache element
457 cache = new CacheElem ;
458
459 // Make list of function projection and normalization integrals
460 for (const auto elm : funcList) {
461 const auto func = static_cast<RooAbsReal*>(elm);
462
463 std::unique_ptr<RooAbsReal> funcInt{func->createIntegral(analVars,rangeName)};
464 if(auto funcRealInt = dynamic_cast<RooRealIntegral*>(funcInt.get())) funcRealInt->setAllowComponentSelection(true);
465 cache->_funcIntList.addOwned(std::move(funcInt));
466 if (normSet && !normSet->empty()) {
467 cache->_funcNormList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(*normSet)});
468 }
469 }
470
471 // Store cache element
472 Int_t code = normIntMgr.setObj(normSet.get(),&analVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
473
474 //cout << "RooRealSumPdf("<<this<<")::getAnalyticalIntegralWN:"<<GetName()<<"("<<allVars<<","<<analVars<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << " -> " << code+1 << endl;
475 return code+1 ;
476}
477
478
479
480
481////////////////////////////////////////////////////////////////////////////////
482/// Implement analytical integrations by deferring integration of component
483/// functions to integrators of components.
484
485double RooRealSumPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* rangeName) const
486{
487 return analyticalIntegralWN(*this, _normIntMgr, _funcList, _coefList, code, normSet2, rangeName, _haveWarned);
488}
489
491 RooArgList const& funcList, RooArgList const& coefList,
492 Int_t code, const RooArgSet* normSet2, const char* rangeName,
493 bool hasWarnedBefore)
494{
495 // Handle trivial passthrough scenario
496 if (code==0) return caller.getVal(normSet2) ;
497
498
499 // WVE needs adaptation for rangeName feature
500 auto* cache = static_cast<CacheElem*>(normIntMgr.getObjByIndex(code-1));
501 if (cache==0) { // revive the (sterilized) cache
502 //cout << "RooRealSumPdf("<<this<<")::analyticalIntegralWN:"<<GetName()<<"("<<code<<","<<(normSet2?*normSet2:RooArgSet())<<","<<(rangeName?rangeName:"<none>") << ": reviving cache "<< endl;
503 RooArgSet vars;
504 caller.getParameters(nullptr, vars);
505 RooArgSet iset = normIntMgr.selectFromSet2(vars, code-1);
506 RooArgSet nset = normIntMgr.selectFromSet1(vars, code-1);
507 RooArgSet dummy;
508 Int_t code2 = caller.getAnalyticalIntegralWN(iset,dummy,&nset,rangeName);
509 R__ASSERT(code==code2); // must have revived the right (sterilized) slot...
510 cache = (CacheElem*) normIntMgr.getObjByIndex(code-1) ;
511 R__ASSERT(cache!=0);
512 }
513
514 double value(0) ;
515
516 // N funcs, N-1 coefficients
517 double lastCoef(1) ;
518 auto funcIt = funcList.begin();
519 auto funcIntIt = cache->_funcIntList.begin();
520 for (const auto coefArg : coefList) {
521 assert(funcIt != funcList.end());
522 const auto coef = static_cast<const RooAbsReal*>(coefArg);
523 const auto func = static_cast<const RooAbsReal*>(*funcIt++);
524 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt++);
525
526 double coefVal = coef->getVal(normSet2) ;
527 if (coefVal) {
528 assert(func);
529 if (normSet2 ==0 || func->isSelectedComp()) {
530 assert(funcInt);
531 value += funcInt->getVal()*coefVal ;
532 }
533 lastCoef -= coef->getVal(normSet2) ;
534 }
535 }
536
537 const bool haveLastCoef = funcList.size() == coefList.size();
538
539 if (!haveLastCoef) {
540 // Add last func with correct coefficient
541 const auto func = static_cast<const RooAbsReal*>(*funcIt);
542 const auto funcInt = static_cast<RooAbsReal*>(*funcIntIt);
543 assert(func);
544
545 if (normSet2 ==0 || func->isSelectedComp()) {
546 assert(funcInt);
547 value += funcInt->getVal()*lastCoef ;
548 }
549
550 // Warn about coefficient degeneration
551 if (!hasWarnedBefore && (lastCoef<0 || lastCoef>1)) {
552 oocoutW(&caller, Eval) << caller.ClassName() << "::evaluate(" << caller.GetName()
553 << " WARNING: sum of FUNC coefficients not in range [0-1], value="
554 << 1-lastCoef << endl ;
555 }
556 }
557
558 double normVal(1) ;
559 if (normSet2 && normSet2->getSize()>0) {
560 normVal = 0 ;
561
562 // N funcs, N-1 coefficients
563 auto funcNormIter = cache->_funcNormList.begin();
564 for (const auto coefAsArg : coefList) {
565 auto coef = static_cast<RooAbsReal*>(coefAsArg);
566 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter++);
567
568 double coefVal = coef->getVal(normSet2);
569 if (coefVal) {
570 assert(funcNorm);
571 normVal += funcNorm->getVal()*coefVal ;
572 }
573 }
574
575 // Add last func with correct coefficient
576 if (!haveLastCoef) {
577 auto funcNorm = static_cast<RooAbsReal*>(*funcNormIter);
578 assert(funcNorm);
579
580 normVal += funcNorm->getVal()*lastCoef;
581 }
582 }
583
584 return value / normVal;
585}
586
587
588////////////////////////////////////////////////////////////////////////////////
589
591{
592 double n = getNorm(nset) ;
593 if (n<0) {
594 logEvalError("Expected number of events is negative") ;
595 }
596 return n ;
597}
598
599
600////////////////////////////////////////////////////////////////////////////////
601
602std::list<double>* RooRealSumPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
603{
604 return binBoundaries(_funcList, obs, xlo, xhi);
605}
606
607
608std::list<double>* RooRealSumPdf::binBoundaries(RooArgList const& funcList, RooAbsRealLValue& obs, double xlo, double xhi)
609{
610 std::list<double>* sumBinB = nullptr;
611 bool needClean(false) ;
612
613 // Loop over components pdf
614 for (auto * func : static_range_cast<RooAbsReal*>(funcList)) {
615
616 list<double>* funcBinB = func->binBoundaries(obs,xlo,xhi) ;
617
618 // Process hint
619 if (funcBinB) {
620 if (!sumBinB) {
621 // If this is the first hint, then just save it
622 sumBinB = funcBinB ;
623 } else {
624
625 std::list<double>* newSumBinB = new list<double>(sumBinB->size()+funcBinB->size()) ;
626
627 // Merge hints into temporary array
628 merge(funcBinB->begin(),funcBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
629
630 // Copy merged array without duplicates to new sumBinBArrau
631 delete sumBinB ;
632 delete funcBinB ;
633 sumBinB = newSumBinB ;
634 needClean = true ;
635 }
636 }
637 }
638
639 // Remove consecutive duplicates
640 if (needClean) {
641 list<double>::iterator new_end = unique(sumBinB->begin(),sumBinB->end()) ;
642 sumBinB->erase(new_end,sumBinB->end()) ;
643 }
644
645 return sumBinB ;
646}
647
648
649
650/// Check if all components that depend on `obs` are binned.
652{
653 return isBinnedDistribution(_funcList, obs);
654}
655
656
658{
659 for (auto* func : static_range_cast<RooAbsReal*>(funcList)) {
660
661 if (func->dependsOn(obs) && !func->isBinnedDistribution(obs)) {
662 return false ;
663 }
664 }
665
666 return true ;
667}
668
669
670
671////////////////////////////////////////////////////////////////////////////////
672
673std::list<double>* RooRealSumPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
674{
675 return plotSamplingHint(_funcList, obs, xlo, xhi);
676}
677
678std::list<double>* RooRealSumPdf::plotSamplingHint(RooArgList const& funcList, RooAbsRealLValue& obs, double xlo, double xhi)
679{
680 std::list<double>* sumHint = nullptr;
681 bool needClean(false) ;
682
683 // Loop over components pdf
684 for (const auto elm : funcList) {
685 auto func = static_cast<RooAbsReal*>(elm);
686
687 list<double>* funcHint = func->plotSamplingHint(obs,xlo,xhi) ;
688
689 // Process hint
690 if (funcHint) {
691 if (!sumHint) {
692
693 // If this is the first hint, then just save it
694 sumHint = funcHint ;
695
696 } else {
697
698 auto* newSumHint = new std::list<double>(sumHint->size()+funcHint->size()) ;
699
700 // the lists must be sorted before merging them
701 funcHint->sort();
702 sumHint->sort();
703 // Merge hints into temporary array
704 merge(funcHint->begin(),funcHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
705
706 // Copy merged array without duplicates to new sumHintArrau
707 delete sumHint ;
708 sumHint = newSumHint ;
709 needClean = true ;
710 }
711 }
712 }
713
714 // Remove consecutive duplicates
715 if (needClean) {
716 sumHint->erase(std::unique(sumHint->begin(),sumHint->end()), sumHint->end()) ;
717 }
718
719 return sumHint ;
720}
721
722
723
724
725////////////////////////////////////////////////////////////////////////////////
726/// Label OK'ed components of a RooRealSumPdf with cache-and-track
727
729{
730 setCacheAndTrackHints(_funcList, trackNodes);
731}
732
733
735{
736 for (const auto sarg : funcList) {
737 if (sarg->canNodeBeCached()==Always) {
738 trackNodes.add(*sarg) ;
739 //cout << "tracking node RealSumPdf component " << sarg->ClassName() << "::" << sarg->GetName() << endl ;
740 }
741 }
742}
743
744
745////////////////////////////////////////////////////////////////////////////////
746/// Customized printing of arguments of a RooRealSumPdf to more intuitively reflect the contents of the
747/// product operator construction
748
749void RooRealSumPdf::printMetaArgs(ostream& os) const
750{
752}
753
754
755void RooRealSumPdf::printMetaArgs(RooArgList const& funcList, RooArgList const& coefList, ostream& os)
756{
757
758 bool first(true) ;
759
760 if (!coefList.empty()) {
761 auto funcIter = funcList.begin();
762
763 for (const auto coef : coefList) {
764 if (!first) {
765 os << " + " ;
766 } else {
767 first = false ;
768 }
769 const auto func = *(funcIter++);
770 os << coef->GetName() << " * " << func->GetName();
771 }
772
773 if (funcIter != funcList.end()) {
774 os << " + [%] * " << (*funcIter)->GetName() ;
775 }
776 } else {
777
778 for (const auto func : funcList) {
779 if (!first) {
780 os << " + " ;
781 } else {
782 first = false ;
783 }
784 os << func->GetName() ;
785 }
786 }
787
788 os << " " ;
789}
790
791std::unique_ptr<RooAbsArg>
793{
794 if (normSet.empty() || selfNormalized()) {
795 return RooAbsPdf::compileForNormSet({}, ctx);
796 }
797 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
798
799 if (ctx.likelihoodMode() && pdfClone->getAttribute("BinnedLikelihood")) {
800
801 // This has to be done before compiling the servers, such that the
802 // RooBinWidthFunctions know to disable themselves.
803 ctx.setBinnedLikelihoodMode(true);
804
805 ctx.markAsCompiled(*pdfClone);
806 ctx.compileServers(*pdfClone, {});
807
808 pdfClone->setAttribute("BinnedLikelihoodActive");
809 // If this is a binned likelihood, we're flagging it in the context.
810 // Then, the RooBinWidthFunctions know that they should not put
811 // themselves in the computation graph. Like this, the pdf values can
812 // directly be interpreted as yields, without multiplying them with the
813 // bin widths again in the NLL. However, the NLL class has to be careful
814 // to only skip the bin with multiplication when there actually were
815 // RooBinWidthFunctions! This is not the case for old workspace before
816 // ROOT 6.26. Therefore, we use the "BinnedLikelihoodActiveYields"
817 // attribute to let the NLL know what it should do.
818 if (ctx.binWidthFuncFlag()) {
819 pdfClone->setAttribute("BinnedLikelihoodActiveYields");
820 }
821 return pdfClone;
822 }
823
824 ctx.compileServers(*pdfClone, {});
825
826 RooArgSet depList;
827 pdfClone->getObservables(&normSet, depList);
828
829 auto newArg = std::make_unique<RooNormalizedPdf>(*pdfClone, depList);
830
831 // The direct servers are this pdf and the normalization integral, which
832 // don't need to be compiled further.
833 for (RooAbsArg *server : newArg->servers()) {
834 ctx.markAsCompiled(*server);
835 }
836 ctx.markAsCompiled(*newArg);
837 newArg->addOwnedComponents(std::move(pdfClone));
838 return newArg;
839}
840
841std::unique_ptr<RooAbsReal> RooRealSumPdf::createExpectedEventsFunc(const RooArgSet *nset) const
842{
843 if (nset == nullptr)
844 return nullptr;
845 return std::unique_ptr<RooAbsReal>{createIntegral(*nset, *getIntegratorConfig(), normRange())};
846}
#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:118
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
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:91
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
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.
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
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
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 ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
bool getForceNumInt() const
Definition RooAbsReal.h:174
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 containing the objects that are both in the cached set 1 with a given index and an i...
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 2 with a given index and an i...
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
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 const & getResult(RooAbsArg const &arg)
Gets the result for the given node using the node name.
std::string getTmpVarName()
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)
Definition DataMap.cxx:22
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...
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.
const RooArgList & funcList() const
~RooRealSumPdf() override
Destructor.
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 computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) 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 void translateImpl(RooFit::Detail::CodeSquashContext &ctx, RooAbsArg const *klass, RooArgList const &funcList, RooArgList const &coefList)
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
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 uint64_t sum(uint64_t i)
Definition Factory.cxx:2345
static void output()