Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAddPdf.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 RooAddPdf
19 \ingroup Roofitcore
20
21RooAddPdf is an efficient implementation of a sum of PDFs of the form
22
23\f[
24 \sum_{i=1}^{n} c_i \cdot \mathrm{PDF}_i
25\f]
26
27or
28\f[
29 c_1\cdot\mathrm{PDF}_1 + c_2\cdot\mathrm{PDF}_2 \; + \; ... \; + \; \left( 1-\sum_{i=1}^{n-1}c_i \right) \cdot \mathrm{PDF}_n
30\f]
31
32The first form is for extended likelihood fits, where the
33expected number of events is \f$ \sum_i c_i \f$. The coefficients \f$ c_i \f$
34can either be explicitly provided, or, if all components support
35extended likelihood fits, they can be calculated from the contribution
36of each PDF to the total expected number of events.
37
38In the second form, the sum of the coefficients is required to be 1 or less,
39and the coefficient of the last PDF is calculated automatically from the condition
40that the sum of all coefficients has to be 1.
41
42### Recursive coefficients
43It is also possible to parameterise the coefficients recursively
44
45\f[
46 \sum_{i=1}^n c_i \prod_{j=1}^{i-1} \left[ (1-c_j) \right] \cdot \mathrm{PDF}_i \\
47 = c_1 \cdot \mathrm{PDF}_1 + (1-c_1)\, c_2 \cdot \mathrm{PDF}_2 + \ldots + (1-c_1)\ldots(1-c_{n-1}) \cdot 1 \cdot \mathrm{PDF}_n \\
48\f]
49
50In this form the sum of the coefficients is always less than 1.0
51for all possible values of the individual coefficients between 0 and 1.
52\note Don't pass the \f$ n^\mathrm{th} \f$ coefficient. It is always 1, since the normalisation condition removes one degree of freedom.
53
54RooAddPdf relies on each component PDF to be normalized and will perform
55no normalization other than calculating the proper last coefficient \f$ c_n \f$, if requested.
56An (enforced) condition for this assumption is that each \f$ \mathrm{PDF}_i \f$ is independent of each \f$ c_i \f$.
57
58## Difference between RooAddPdf / RooRealSumFunc / RooRealSumPdf
59- RooAddPdf is a PDF of PDFs, *i.e.* its components need to be normalised and non-negative.
60- RooRealSumPdf is a PDF of functions, *i.e.*, its components can be negative, but their sum cannot be. The normalisation
61 is computed automatically, unless the PDF is extended (see above).
62- RooRealSumFunc is a sum of functions. It is neither normalised, nor need it be positive.
63
64*/
65
66#include "RooAddPdf.h"
67
68#include "RooAddGenContext.h"
69#include "RooBatchCompute.h"
70#include "RooDataSet.h"
71#include "RooGlobalFunc.h"
72#include "RooNaNPacker.h"
73#include "RooRealProxy.h"
74#include "RooRealVar.h"
75#include "RooRealConstant.h"
76#include "RooRealIntegral.h"
78
79#include <algorithm>
80#include <memory>
81#include <sstream>
82#include <set>
83
84using namespace std;
85
87
88
89////////////////////////////////////////////////////////////////////////////////
90/// Dummy constructor
91
92RooAddPdf::RooAddPdf(const char *name, const char *title) :
93 RooAbsPdf(name,title),
94 _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,false,false),
95 _projCacheMgr(this,10),
96 _pdfList("!pdfs","List of PDFs",this),
97 _coefList("!coefficients","List of coefficients",this),
98 _coefErrCount{_errorCount}
99{
101}
102
103
105
106 // Two pdfs with the same name are only allowed in the input list if they are
107 // actually the same object.
108 using PdfInfo = std::pair<std::string,RooAbsArg*>;
109 std::set<PdfInfo> seen;
110 for(auto const& pdf : _pdfList) {
111 PdfInfo elem{pdf->GetName(), pdf};
112 auto comp = [&](PdfInfo const& p){ return p.first == elem.first && p.second != elem.second; };
113 auto found = std::find_if(seen.begin(), seen.end(), comp);
114 if(found != seen.end()) {
115 std::stringstream errorMsg;
116 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
117 << ") pdf list contains pdfs with duplicate name \"" << pdf->GetName() << "\"."
118 << std::endl;
119 coutE(InputArguments) << errorMsg.str();
120 throw std::invalid_argument(errorMsg.str().c_str());
121 }
122 seen.insert(elem);
123 }
124
125 _coefCache.resize(_pdfList.size());
126}
127
128
129////////////////////////////////////////////////////////////////////////////////
130/// Constructor with two PDFs and one coefficient
131
132RooAddPdf::RooAddPdf(const char *name, const char *title,
133 RooAbsPdf& pdf1, RooAbsPdf& pdf2, RooAbsReal& coef1) :
134 RooAddPdf(name, title)
135{
136 _pdfList.add(pdf1) ;
137 _pdfList.add(pdf2) ;
138 _coefList.add(coef1) ;
139
141}
142
143
144////////////////////////////////////////////////////////////////////////////////
145/// Generic constructor from list of PDFs and list of coefficients.
146/// Each pdf list element (i) is paired with coefficient list element (i).
147/// The number of coefficients must be either equal to the number of PDFs,
148/// in which case extended MLL fitting is enabled, or be one less.
149///
150/// All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal
151///
152/// If the recursiveFraction flag is true, the coefficients are interpreted as recursive
153/// coefficients as explained in the class description.
154
155RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, bool recursiveFractions) :
156 RooAddPdf(name,title)
157{
158 _recursive = recursiveFractions;
159
160 if (inPdfList.size()>inCoefList.size()+1 || inPdfList.size()<inCoefList.size()) {
161 std::stringstream errorMsg;
162 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
163 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1." << endl ;
164 coutE(InputArguments) << errorMsg.str();
165 throw std::invalid_argument(errorMsg.str().c_str());
166 }
167
168 if (recursiveFractions && inPdfList.size()!=inCoefList.size()+1) {
169 std::stringstream errorMsg;
170 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
171 << "): Recursive fractions option can only be used if Npdf=Ncoef+1." << endl;
172 coutE(InputArguments) << errorMsg.str();
173 throw std::invalid_argument(errorMsg.str());
174 }
175
176 // Constructor with N PDFs and N or N-1 coefs
177 RooArgList partinCoefList ;
178
179 bool first(true) ;
180
181 for (auto i = 0u; i < inCoefList.size(); ++i) {
182 auto coef = dynamic_cast<RooAbsReal*>(inCoefList.at(i));
183 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(i));
184 if (inPdfList.at(i) == nullptr) {
185 std::stringstream errorMsg;
186 errorMsg << "RooAddPdf::RooAddPdf(" << GetName()
187 << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
188 coutE(InputArguments) << errorMsg.str();
189 throw std::invalid_argument(errorMsg.str());
190 }
191 if (!coef) {
192 std::stringstream errorMsg;
193 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") coefficient " << (coef ? coef->GetName() : "") << " is not of type RooAbsReal, ignored" << endl ;
194 coutE(InputArguments) << errorMsg.str();
195 throw std::invalid_argument(errorMsg.str());
196 }
197 if (!pdf) {
198 std::stringstream errorMsg;
199 errorMsg << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored" << endl ;
200 coutE(InputArguments) << errorMsg.str();
201 throw std::invalid_argument(errorMsg.str());
202 }
203 _pdfList.add(*pdf) ;
204
205 // Process recursive fraction mode separately
206 if (recursiveFractions) {
207 partinCoefList.add(*coef) ;
208 if (first) {
209
210 // The first fraction is the first plain fraction
211 first = false ;
212 _coefList.add(*coef) ;
213
214 } else {
215
216 // The i-th recursive fraction = (1-f1)*(1-f2)*...(fi) and is calculated from the list (f1,...,fi) by RooRecursiveFraction)
217 RooAbsReal* rfrac = new RooRecursiveFraction(Form("%s_recursive_fraction_%s",GetName(),pdf->GetName()),"Recursive Fraction",partinCoefList) ;
218 addOwnedComponents(*rfrac) ;
219 _coefList.add(*rfrac) ;
220
221 }
222
223 } else {
224 _coefList.add(*coef) ;
225 }
226 }
227
228 if (inPdfList.size() == inCoefList.size() + 1) {
229 auto pdf = dynamic_cast<RooAbsPdf*>(inPdfList.at(inCoefList.size()));
230
231 if (!pdf) {
232 coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") last argument " << inPdfList.at(inCoefList.size())->GetName() << " is not of type RooAbsPdf." << endl ;
233 throw std::invalid_argument("Last argument for RooAddPdf is not a PDF.");
234 }
235 _pdfList.add(*pdf) ;
236
237 // Process recursive fractions mode. Above, we verified that we don't have a last coefficient
238 if (recursiveFractions) {
239
240 // The last recursive fraction = (1-f1)*(1-f2)*...(1-fN) and is calculated from the list (f1,...,fN,1) by RooRecursiveFraction
241 partinCoefList.add(RooFit::RooConst(1)) ;
242 RooAbsReal* rfrac = new RooRecursiveFraction(Form("%s_recursive_fraction_%s",GetName(),pdf->GetName()),"Recursive Fraction",partinCoefList) ;
243 addOwnedComponents(*rfrac) ;
244 _coefList.add(*rfrac) ;
245
246 // In recursive mode we always have Ncoef=Npdf, since we added it just above
247 _haveLastCoef=true ;
248 }
249
250 } else {
251 _haveLastCoef=true ;
252 }
253
255}
256
257
258////////////////////////////////////////////////////////////////////////////////
259/// Generic constructor from list of extended PDFs. There are no coefficients as the expected
260/// number of events from each components determine the relative weight of the PDFs.
261///
262/// All PDFs must inherit from RooAbsPdf.
263
264RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList& inPdfList) :
265 RooAddPdf(name,title)
266{
267 _allExtendable = true;
268
269 // Constructor with N PDFs
270 for (const auto pdfArg : inPdfList) {
271 auto pdf = dynamic_cast<const RooAbsPdf*>(pdfArg);
272
273 if (!pdf) {
274 coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << (pdf ? pdf->GetName() : "") << " is not of type RooAbsPdf, ignored" << endl ;
275 continue ;
276 }
277 if (!pdf->canBeExtended()) {
278 coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName() << " is not extendable, ignored" << endl ;
279 continue ;
280 }
281 _pdfList.add(*pdf) ;
282 }
283
285}
286
287
288////////////////////////////////////////////////////////////////////////////////
289/// Copy constructor
290
291RooAddPdf::RooAddPdf(const RooAddPdf& other, const char* name) :
292 RooAbsPdf(other,name),
293 _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
294 _refCoefRangeName((TNamed*)other._refCoefRangeName),
295 _projectCoefs(other._projectCoefs),
296 _projCacheMgr(other._projCacheMgr,this),
297 _codeReg(other._codeReg),
298 _pdfList("!pdfs",this,other._pdfList),
299 _coefList("!coefficients",this,other._coefList),
300 _haveLastCoef(other._haveLastCoef),
301 _allExtendable(other._allExtendable),
302 _recursive(other._recursive)
303{
307}
308
309
310////////////////////////////////////////////////////////////////////////////////
311/// By default the interpretation of the fraction coefficients is
312/// performed in the contextual choice of observables. This makes the
313/// shape of the p.d.f explicitly dependent on the choice of
314/// observables. This method instructs RooAddPdf to freeze the
315/// interpretation of the coefficients to be done in the given set of
316/// observables. If frozen, fractions are automatically transformed
317/// from the reference normalization set to the contextual normalization
318/// set by ratios of integrals.
319
321{
322 if (refCoefNorm.empty()) {
323 _projectCoefs = false ;
324 return ;
325 }
326 _projectCoefs = true ;
327
329 _refCoefNorm.add(refCoefNorm) ;
330
332}
333
334
335
336////////////////////////////////////////////////////////////////////////////////
337/// By default, fraction coefficients are assumed to refer to the default
338/// fit range. This makes the shape of a RooAddPdf
339/// explicitly dependent on the range of the observables. Calling this function
340/// allows for a range-independent definition of the fractions, because it
341/// ties all coefficients to the given
342/// named range. If the normalisation range is different
343/// from this reference range, the appropriate fraction coefficients
344/// are automatically calculated from the reference fractions by
345/// integrating over the ranges, and comparing these integrals.
346
347void RooAddPdf::fixCoefRange(const char* rangeName)
348{
350 if (_refCoefRangeName) _projectCoefs = true ;
351}
352
353
354
355////////////////////////////////////////////////////////////////////////////////
356/// Retrieve cache element for the computation of the PDF normalisation.
357/// \param[in] nset Current normalisation set (integration over these variables yields 1).
358/// \param[in] iset Integration set. Variables to be integrated over (if integrations are performed).
359/// \param[in] rangeName Reference range for the integrals.
360///
361/// If a cache element does not exist, create and fill it on the fly. The cache also contains
362/// - Supplemental normalization terms (in case not all added p.d.f.s have the same observables)
363/// - Projection integrals to calculate transformed fraction coefficients when a frozen reference frame is provided
364/// - Projection integrals for similar transformations when a frozen reference range is provided.
365
366RooAddPdf::CacheElem* RooAddPdf::getProjCache(const RooArgSet* nset, const RooArgSet* iset, const char* rangeName) const
367{
368
369 // Check if cache already exists
370 CacheElem* cache = (CacheElem*) _projCacheMgr.getObj(nset,iset,0,rangeName) ;
371 if (cache) {
372 return cache ;
373 }
374
375 //Create new cache
376 cache = new CacheElem ;
377
378 // *** PART 1 : Create supplemental normalization list ***
379
380 // Retrieve the combined set of dependents of this PDF ;
381 auto fullDepList = std::unique_ptr<RooArgSet>{getObservables(nset)} ;
382 if (iset) {
383 fullDepList->remove(*iset,true,true) ;
384 }
385
386 // Fill with dummy unit RRVs for now
387 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
388 auto pdf = static_cast<const RooAbsPdf *>(_pdfList.at(i));
389 auto coef = static_cast<const RooAbsReal*>(_coefList.at(i));
390
391 // Start with full list of dependents
392 RooArgSet supNSet(*fullDepList) ;
393
394 // Remove PDF dependents
395 if (auto pdfDeps = std::unique_ptr<RooArgSet>{pdf->getObservables(nset)}) {
396 supNSet.remove(*pdfDeps,true,true) ;
397 }
398
399 // Remove coef dependents
400 if (auto coefDeps = std::unique_ptr<RooArgSet>{coef ? coef->getObservables(nset) : nullptr}) {
401 supNSet.remove(*coefDeps,true,true) ;
402 }
403
404 RooAbsReal* snorm ;
405 auto name = std::string(GetName()) + "_" + pdf->GetName() + "_SupNorm";
406 cache->_needSupNorm = false ;
407 if (!supNSet.empty()) {
408 snorm = new RooRealIntegral(name.c_str(),"Supplemental normalization integral",RooRealConstant::value(1.0),supNSet) ;
409 cxcoutD(Caching) << "RooAddPdf " << GetName() << " making supplemental normalization set " << supNSet << " for pdf component " << pdf->GetName() << endl ;
410 cache->_needSupNorm = true ;
411 } else {
412 snorm = new RooRealVar(name.c_str(),"Unit Supplemental normalization integral",1.0) ;
413 }
414 cache->_suppNormList.addOwned(*snorm) ;
415 }
416
417 if (_verboseEval>1) {
418 cxcoutD(Caching) << "RooAddPdf::syncSuppNormList(" << GetName() << ") synching supplemental normalization list for norm" << (nset?*nset:RooArgSet()) << endl ;
419 if dologD(Caching) {
420 cache->_suppNormList.Print("v") ;
421 }
422 }
423
424
425 // *** PART 2 : Create projection coefficients ***
426
427// cout << " this = " << this << " (" << GetName() << ")" << endl ;
428// cout << "projectCoefs = " << (_projectCoefs?"T":"F") << endl ;
429// cout << "_normRange.Length() = " << _normRange.Length() << endl ;
430
431 // If no projections required stop here
432 if (!_projectCoefs && !rangeName) {
433 _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
434// cout << " no projection required" << endl ;
435 return cache ;
436 }
437
438
439// cout << "calculating projection" << endl ;
440
441 // Reduce iset/nset to actual dependents of this PDF
442 auto nset2 = std::unique_ptr<RooArgSet>{nset ? getObservables(nset) : new RooArgSet()} ;
443 cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset = " << (nset?*nset:RooArgSet()) << " nset2 = " << *nset2 << endl ;
444
445 if (nset2->empty() && !_refCoefNorm.empty()) {
446 //cout << "WVE: evaluating RooAddPdf without normalization, but have reference normalization for coefficient definition" << endl ;
447
448 nset2->add(_refCoefNorm) ;
449 if (_refCoefRangeName) {
451 }
452 }
453
454
455 // Check if requested transformation is not identity
456 if (!nset2->equals(_refCoefNorm) || _refCoefRangeName !=0 || rangeName !=0 || _normRange.Length()>0) {
457
458 cxcoutD(Caching) << "ALEX: RooAddPdf::syncCoefProjList(" << GetName() << ") projecting coefficients from "
459 << *nset2 << (rangeName?":":"") << (rangeName?rangeName:"")
460 << " to " << ((!_refCoefNorm.empty())?_refCoefNorm:*nset2) << (_refCoefRangeName?":":"") << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"") << endl ;
461
462 // Recalculate projection integrals of PDFs
463 for (auto arg : _pdfList) {
464 auto thePdf = static_cast<const RooAbsPdf*>(arg);
465
466 // Calculate projection integral
467 RooAbsReal* pdfProj ;
468 if (!nset2->equals(_refCoefNorm)) {
469 pdfProj = thePdf->createIntegral(*nset2,_refCoefNorm,_normRange.Length()>0?_normRange.Data():0) ;
470 pdfProj->setOperMode(operMode()) ;
471 cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset2(" << *nset2 << ")!=_refCoefNorm(" << _refCoefNorm << ") --> pdfProj = " << pdfProj->GetName() << endl ;
472 } else {
473 auto name = std::string(GetName()) + "_" + thePdf->GetName() + "_ProjectNorm";
474 pdfProj = new RooRealVar(name.c_str(),"Unit Projection normalization integral",1.0) ;
475 cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset2(" << *nset2 << ")==_refCoefNorm("
476 << _refCoefNorm << ") --> pdfProj = " << pdfProj->GetName() << std::endl;
477 }
478
479 cache->_projList.addOwned(*pdfProj) ;
480 cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") PP = " << pdfProj->GetName() << endl ;
481
482 // Calculation optional supplemental normalization term
483 RooArgSet supNormSet(_refCoefNorm) ;
484 auto deps = std::unique_ptr<RooArgSet>{thePdf->getParameters(RooArgSet())} ;
485 supNormSet.remove(*deps,true,true) ;
486
487 RooAbsReal* snorm ;
488 auto name = std::string(GetName()) + "_" + thePdf->GetName() + "_ProjSupNorm";
489 if (!supNormSet.empty() && !nset2->equals(_refCoefNorm) ) {
490 snorm = new RooRealIntegral(name.c_str(),"Projection Supplemental normalization integral",
491 RooRealConstant::value(1.0),supNormSet) ;
492 } else {
493 snorm = new RooRealVar(name.c_str(),"Unit Projection Supplemental normalization integral",1.0) ;
494 }
495 cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") SN = " << snorm->GetName() << endl ;
496 cache->_suppProjList.addOwned(*snorm) ;
497
498 // Calculate reference range adjusted projection integral
499 RooAbsReal* rangeProj1 ;
500
501 // cout << "ALEX >>>> RooAddPdf(" << GetName() << ")::getPC _refCoefRangeName WVE = "
502// <<(_refCoefRangeName?":":"") << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"")
503// <<" _refCoefRangeName AK = " << (_refCoefRangeName?_refCoefRangeName->GetName():"")
504// << " && _refCoefNorm" << _refCoefNorm << " with size = _refCoefNorm.size() " << _refCoefNorm.size() << endl ;
505
506 // Check if _refCoefRangeName is identical to default range for all observables,
507 // If so, substitute by unit integral
508
509 // ----------
510 auto tmpObs = std::unique_ptr<RooArgSet>{thePdf->getObservables(_refCoefNorm)} ;
511 bool allIdent = true ;
512 for (auto const& obsArg : *tmpObs) {
513 RooRealVar* rvarg = dynamic_cast<RooRealVar*>(obsArg) ;
514 if (rvarg) {
515 if (rvarg->getMin(RooNameReg::str(_refCoefRangeName))!=rvarg->getMin() ||
516 rvarg->getMax(RooNameReg::str(_refCoefRangeName))!=rvarg->getMax()) {
517 allIdent=false ;
518 }
519 }
520 }
521 // -------------
522
523 if (_refCoefRangeName && !_refCoefNorm.empty() && !allIdent) {
524
525
526 auto tmp = std::unique_ptr<RooArgSet>{thePdf->getObservables(_refCoefNorm)} ;
527 rangeProj1 = thePdf->createIntegral(*tmp,*tmp,RooNameReg::str(_refCoefRangeName)) ;
528
529 //rangeProj1->setOperMode(operMode()) ;
530
531 } else {
532
533 auto theName = std::string(GetName()) + "_" + thePdf->GetName() + "_RangeNorm1";
534 rangeProj1 = new RooRealVar(theName.c_str(),"Unit range normalization integral",1.0) ;
535
536 }
537 cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") R1 = " << rangeProj1->GetName() << endl ;
538 cache->_refRangeProjList.addOwned(*rangeProj1) ;
539
540
541 // Calculate range adjusted projection integral
542 RooAbsReal* rangeProj2 ;
543 cxcoutD(Caching) << "RooAddPdf::syncCoefProjList(" << GetName() << ") rangename = " << (rangeName?rangeName:"<null>")
544 << " nset = " << (nset?*nset:RooArgSet()) << endl ;
545 if (rangeName && !_refCoefNorm.empty()) {
546
547 rangeProj2 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,rangeName) ;
548 //rangeProj2->setOperMode(operMode()) ;
549
550 } else if (_normRange.Length()>0) {
551
552 auto tmp = std::unique_ptr<RooArgSet>{thePdf->getObservables(_refCoefNorm)} ;
553 rangeProj2 = thePdf->createIntegral(*tmp,*tmp,_normRange.Data()) ;
554
555 } else {
556
557 auto theName = std::string(GetName()) + "_" + thePdf->GetName() + "_RangeNorm2";
558 rangeProj2 = new RooRealVar(theName.c_str(),"Unit range normalization integral",1.0) ;
559
560 }
561 cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") R2 = " << rangeProj2->GetName() << endl ;
562 cache->_rangeProjList.addOwned(*rangeProj2) ;
563
564 }
565
566 }
567
568 _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
569
570 return cache ;
571}
572
573
574////////////////////////////////////////////////////////////////////////////////
575/// Update the coefficient values in the given cache element: calculate new remainder
576/// fraction, normalize fractions obtained from extended ML terms to unity, and
577/// multiply the various range and dimensional corrections needed in the
578/// current use context.
579
580void RooAddPdf::updateCoefficients(CacheElem& cache, const RooArgSet* nset) const
581{
582 // Since this function updates the cache, it obviously needs write access:
583 auto& myCoefCache = const_cast<std::vector<double>&>(_coefCache);
584 myCoefCache.resize(_haveLastCoef ? _coefList.size() : _pdfList.size(), 0.);
585
586 // Straight coefficients
587 if (_allExtendable) {
588
589 // coef[i] = expectedEvents[i] / SUM(expectedEvents)
590 double coefSum(0) ;
591 std::size_t i = 0;
592 for (auto arg : _pdfList) {
593 auto pdf = static_cast<RooAbsPdf*>(arg);
594 myCoefCache[i] = pdf->expectedEvents(!_refCoefNorm.empty()?&_refCoefNorm:nset) ;
595 coefSum += myCoefCache[i] ;
596 i++ ;
597 }
598
599 if (coefSum==0.) {
600 coutW(Eval) << "RooAddPdf::updateCoefCache(" << GetName() << ") WARNING: total number of expected events is 0" << endl ;
601 } else {
602 for (std::size_t j=0; j < _pdfList.size(); j++) {
603 myCoefCache[j] /= coefSum ;
604 }
605 }
606
607 } else {
608 if (_haveLastCoef) {
609
610 // coef[i] = coef[i] / SUM(coef)
611 double coefSum(0) ;
612 std::size_t i=0;
613 for (auto coefArg : _coefList) {
614 auto coef = static_cast<RooAbsReal*>(coefArg);
615 myCoefCache[i] = coef->getVal(nset) ;
616 coefSum += myCoefCache[i++];
617 }
618 if (coefSum==0.) {
619 coutW(Eval) << "RooAddPdf::updateCoefCache(" << GetName() << ") WARNING: sum of coefficients is zero 0" << endl ;
620 } else {
621 const double invCoefSum = 1./coefSum;
622 for (std::size_t j=0; j < _coefList.size(); j++) {
623 myCoefCache[j] *= invCoefSum;
624 }
625 }
626 } else {
627
628 // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
629 double lastCoef(1) ;
630 std::size_t i=0;
631 for (auto coefArg : _coefList) {
632 auto coef = static_cast<RooAbsReal*>(coefArg);
633 myCoefCache[i] = coef->getVal(nset) ;
634 lastCoef -= myCoefCache[i++];
635 }
636 myCoefCache[_coefList.size()] = lastCoef ;
637
638 // Treat coefficient degeneration
639 const float coefDegen = lastCoef < 0. ? -lastCoef : (lastCoef > 1. ? lastCoef - 1. : 0.);
640 if (coefDegen > 1.E-5) {
641 myCoefCache[_coefList.size()] = RooNaNPacker::packFloatIntoNaN(100.f*coefDegen);
642
643 std::stringstream msg;
644 if (_coefErrCount-->0) {
645 msg << "RooAddPdf::updateCoefCache(" << GetName()
646 << " WARNING: sum of PDF coefficients not in range [0-1], value="
647 << 1-lastCoef ;
648 if (_coefErrCount==0) {
649 msg << " (no more will be printed)" ;
650 }
651 coutW(Eval) << msg.str() << std::endl;
652 }
653 }
654 }
655 }
656
657
658 // Stop here if not projection is required or needed
659 if ((!_projectCoefs && _normRange.Length()==0) || cache._projList.empty()) {
660 return ;
661 }
662
663 // Adjust coefficients for given projection
664 double coefSum(0) ;
665 {
667
668 for (std::size_t i = 0; i < _pdfList.size(); i++) {
669
670 RooAbsReal* pp = ((RooAbsReal*)cache._projList.at(i)) ;
671 RooAbsReal* sn = ((RooAbsReal*)cache._suppProjList.at(i)) ;
672 RooAbsReal* r1 = ((RooAbsReal*)cache._refRangeProjList.at(i)) ;
673 RooAbsReal* r2 = ((RooAbsReal*)cache._rangeProjList.at(i)) ;
674
675 double proj = pp->getVal()/sn->getVal()*(r2->getVal()/r1->getVal()) ;
676
677 myCoefCache[i] *= proj ;
678 coefSum += myCoefCache[i] ;
679 }
680 }
681
682
684 for (std::size_t i=0; i < _pdfList.size(); ++i) {
685 ccoutD(Caching) << " ALEX: POST-SYNC coef[" << i << "] = " << myCoefCache[i]
686 << " ( _coefCache[i]/coefSum = " << myCoefCache[i]*coefSum << "/" << coefSum << " ) "<< endl ;
687 }
688 }
689
690 if (coefSum==0.) {
691 coutE(Eval) << "RooAddPdf::updateCoefCache(" << GetName() << ") sum of coefficients is zero." << endl ;
692 }
693
694 for (std::size_t i=0; i < _pdfList.size(); i++) {
695 myCoefCache[i] /= coefSum ;
696 }
697
698}
699
700////////////////////////////////////////////////////////////////////////////////
701/// Look up projection cache and per-PDF norm sets. If a PDF doesn't have a special
702/// norm set, use the `defaultNorm`. If `defaultNorm == nullptr`, use the member
703/// _normSet.
704std::pair<const RooArgSet*, RooAddPdf::CacheElem*> RooAddPdf::getNormAndCache(const RooArgSet* nset) const {
705
706 // Treat empty normalization set and nullptr the same way.
707 if(nset && nset->empty()) nset = nullptr;
708
709 if (nset == nullptr) {
710 if (!_refCoefNorm.empty()) {
711 nset = &_refCoefNorm ;
712 }
713 }
714
715 // A RooAddPdf needs to have a normalization set defined, otherwise its
716 // coefficient will not be uniquely defined. Its shape depends on the
717 // normalization provided. Un-normalized calls to RooAddPdf can happen in
718 // Roofit, when printing the pdf's or when computing integrals. In these case,
719 // if the pdf has a normalization set previously defined (i.e. stored as a
720 // datamember in _copyOfLastNormSet) it should use it by default when the pdf
721 // is evaluated without passing a normalizations set (in pdf->getVal(nullptr) )
722 // In the case of no pre-defined normalization set exists, a warning will be
723 // produced, since the obtained value will be arbitrary. Note that to avoid
724 // unnecessary warning messages, when calling RooAbsPdf::printValue or
725 // RooAbsPdf::graphVizTree, the printing of the warning messages for the
726 // RooFit::Eval topic is explicitly disabled.
727 {
728 // If nset is still nullptr, get the pointer to a copy of the last-used
729 // normalization set. It nset is not nullptr, check whether the copy of
730 // the last-used normalization set needs an update.
731 if(nset == nullptr) {
732 nset = _copyOfLastNormSet.get();
733 } else if(nset->uniqueId() != _idOfLastUsedNormSet) {
734 _copyOfLastNormSet = std::make_unique<const RooArgSet>(*nset);
736 }
737
738 // If nset is STILL nullptr, print a warning.
739 if (nset == nullptr) {
740 oocoutW(this, Eval) << "Evaluating RooAddPdf without a defined normalization set. This can lead to ambiguos "
741 "coefficients definition and incorrect results."
742 << " Use RooAddPdf::fixCoefNormalization(nset) to provide a normalization set for "
743 "defining uniquely RooAddPdf coefficients!"
744 << std::endl;
745 }
746 }
747
748
749 CacheElem* cache = getProjCache(nset) ;
750 updateCoefficients(*cache,nset) ;
751
752 return {nset, cache};
753}
754
755
756////////////////////////////////////////////////////////////////////////////////
757/// Calculate and return the current value
758
759double RooAddPdf::getValV(const RooArgSet* normSet) const
760{
761 auto normAndCache = getNormAndCache(normSet);
762 const RooArgSet* nset = normAndCache.first;
763 CacheElem* cache = normAndCache.second;
764
765 // Process change in last data set used
766 bool nsetChanged(false) ;
767 if (nset!=_normSet || _norm==0) {
768 nsetChanged = syncNormalization(nset) ;
769 }
770
771 // Do running sum of coef/pdf pairs, calculate lastCoef.
772 if (isValueDirty() || nsetChanged) {
773 _value = 0.0;
774
775 for (unsigned int i=0; i < _pdfList.size(); ++i) {
776 const auto& pdf = static_cast<RooAbsPdf&>(_pdfList[i]);
777 double snormVal = 1.;
778 if (cache->_needSupNorm) {
779 snormVal = ((RooAbsReal*)cache->_suppNormList.at(i))->getVal();
780 }
781
782 double pdfVal = pdf.getVal(nset);
783 if (pdf.isSelectedComp()) {
784 _value += pdfVal*_coefCache[i]/snormVal;
785 }
786 }
788 }
789
790 return _value;
791}
792
793
794////////////////////////////////////////////////////////////////////////////////
795/// Compute addition of PDFs in batches.
796void RooAddPdf::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
797{
800 CacheElem* cache = getNormAndCache(nullptr).second;
801 for (unsigned int pdfNo = 0; pdfNo < _pdfList.size(); ++pdfNo)
802 {
803 auto pdf = static_cast<RooAbsPdf*>(&_pdfList[pdfNo]);
804 if (pdf->isSelectedComp())
805 {
806 pdfs.push_back(dataMap.at(pdf));
807 coefs.push_back(_coefCache[pdfNo] / (cache->_needSupNorm ?
808 static_cast<RooAbsReal*>(cache->_suppNormList.at(pdfNo))->getVal() : 1) );
809 }
810 }
812 dispatch->compute(stream, RooBatchCompute::AddPdf, output, nEvents, pdfs, coefs);
813}
814
815
816////////////////////////////////////////////////////////////////////////////////
817/// Reset error counter to given value, limiting the number
818/// of future error messages for this pdf to 'resetValue'
819
821{
823 _coefErrCount = resetValue ;
824}
825
826
827
828////////////////////////////////////////////////////////////////////////////////
829/// Check if PDF is valid for given normalization set.
830/// Coeffient and PDF must be non-overlapping, but pdf-coefficient
831/// pairs may overlap each other
832
834{
835 bool ret(false) ;
836
837 // There may be fewer coefficients than PDFs.
838 std::size_t end = std::min(_pdfList.size(), _coefList.size());
839 for (std::size_t i = 0; i < end; ++i) {
840 auto pdf = static_cast<const RooAbsPdf *>(_pdfList.at(i));
841 auto coef = static_cast<const RooAbsReal*>(_coefList.at(i));
842 if (pdf->observableOverlaps(nset,*coef)) {
843 coutE(InputArguments) << "RooAddPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName()
844 << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
845 ret = true ;
846 }
847 }
848
849 return ret ;
850}
851
852
853////////////////////////////////////////////////////////////////////////////////
854/// Determine which part (if any) of given integral can be performed analytically.
855/// If any analytical integration is possible, return integration scenario code
856///
857/// RooAddPdf queries each component PDF for its analytical integration capability of the requested
858/// set ('allVars'). It finds the largest common set of variables that can be integrated
859/// by all components. If such a set exists, it reconfirms that each component is capable of
860/// analytically integrating the common set, and combines the components individual integration
861/// codes into a single integration code valid for RooAddPdf.
862
864 const RooArgSet* normSet, const char* rangeName) const
865{
866
867 RooArgSet allAnalVars(*std::unique_ptr<RooArgSet>{getObservables(allVars)}) ;
868
869 Int_t n(0) ;
870
871 // First iteration, determine what each component can integrate analytically
872 for (const auto pdfArg : _pdfList) {
873 auto pdf = static_cast<const RooAbsPdf *>(pdfArg);
874 RooArgSet subAnalVars ;
875 pdf->getAnalyticalIntegralWN(allVars,subAnalVars,normSet,rangeName) ;
876
877 // Observables that cannot be integrated analytically by this component are dropped from the common list
878 for (const auto arg : allVars) {
879 if (!subAnalVars.find(arg->GetName()) && pdf->dependsOn(*arg)) {
880 allAnalVars.remove(*arg,true,true) ;
881 }
882 }
883 n++ ;
884 }
885
886 // If no observables can be integrated analytically, return code 0 here
887 if (allAnalVars.empty()) {
888 return 0 ;
889 }
890
891
892 // Now retrieve codes for integration over common set of analytically integrable observables for each component
893 n=0 ;
894 std::vector<Int_t> subCode(_pdfList.size());
895 bool allOK(true) ;
896 for (const auto arg : _pdfList) {
897 auto pdf = static_cast<const RooAbsPdf *>(arg);
898 RooArgSet subAnalVars ;
899 auto allAnalVars2 = std::unique_ptr<RooArgSet>{pdf->getObservables(allAnalVars)} ;
900 subCode[n] = pdf->getAnalyticalIntegralWN(*allAnalVars2,subAnalVars,normSet,rangeName) ;
901 if (subCode[n]==0 && !allAnalVars2->empty()) {
902 coutE(InputArguments) << "RooAddPdf::getAnalyticalIntegral(" << GetName() << ") WARNING: component PDF " << pdf->GetName()
903 << " advertises inconsistent set of integrals (e.g. (X,Y) but not X or Y individually."
904 << " Distributed analytical integration disabled. Please fix PDF" << endl ;
905 allOK = false ;
906 }
907 n++ ;
908 }
909 if (!allOK) {
910 return 0 ;
911 }
912
913 // Mare all analytically integrated observables as such
914 analVars.add(allAnalVars) ;
915
916 // Store set of variables analytically integrated
917 RooArgSet* intSet = new RooArgSet(allAnalVars) ;
918 Int_t masterCode = _codeReg.store(subCode,intSet)+1 ;
919
920 return masterCode ;
921}
922
923
924
925////////////////////////////////////////////////////////////////////////////////
926/// Return analytical integral defined by given scenario code
927
928Double_t RooAddPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
929{
930 // WVE needs adaptation to handle new rangeName feature
931 if (code==0) {
932 return getVal(normSet) ;
933 }
934
935 // Retrieve analytical integration subCodes and set of observabels integrated over
936 RooArgSet* intSet ;
937 const std::vector<Int_t>& subCode = _codeReg.retrieve(code-1,intSet) ;
938 if (subCode.empty()) {
939 coutE(InputArguments) << "RooAddPdf::analyticalIntegral(" << GetName() << "): ERROR unrecognized integration code, " << code << endl ;
940 assert(0) ;
941 }
942
943 cxcoutD(Caching) << "RooAddPdf::aiWN(" << GetName() << ") calling getProjCache with nset = " << (normSet?*normSet:RooArgSet()) << endl ;
944
945 if ((normSet==0 || normSet->empty()) && !_refCoefNorm.empty()) {
946// cout << "WVE integration of RooAddPdf without normalization, but have reference set, using ref set for normalization" << endl ;
947 normSet = &_refCoefNorm ;
948 }
949
950 CacheElem* cache = getProjCache(normSet,intSet,0) ; // WVE rangename here?
951 updateCoefficients(*cache,normSet) ;
952
953 // Calculate the current value of this object
954 double value(0) ;
955
956 // Do running sum of coef/pdf pairs, calculate lastCoef.
957 double snormVal ;
958
959 //cout << "ROP::aIWN updateCoefCache with rangeName = " << (rangeName?rangeName:"<null>") << endl ;
960 RooArgList* snormSet = (!cache->_suppNormList.empty()) ? &cache->_suppNormList : nullptr ;
961 for (std::size_t i = 0; i < _pdfList.size(); ++i ) {
962 auto pdf = static_cast<const RooAbsPdf*>(_pdfList.at(i));
963
964 if (_coefCache[i]) {
965 snormVal = snormSet ? ((RooAbsReal*) cache->_suppNormList.at(i))->getVal() : 1.0 ;
966
967 // WVE swap this?
968 double val = pdf->analyticalIntegralWN(subCode[i],normSet,rangeName) ;
969 if (pdf->isSelectedComp()) {
970 value += val*_coefCache[i]/snormVal ;
971 }
972 }
973 }
974
975 return value ;
976}
977
978
979
980////////////////////////////////////////////////////////////////////////////////
981/// Return the number of expected events, which is either the sum of all coefficients
982/// or the sum of the components extended terms, multiplied with the fraction that
983/// is in the current range w.r.t the reference range
984
986{
987 double expectedTotal{0.0};
988
989 cxcoutD(Caching) << "RooAddPdf::expectedEvents(" << GetName() << ") calling getProjCache with nset = " << (nset?*nset:RooArgSet()) << endl ;
990 CacheElem& cache = *getProjCache(nset) ;
991 updateCoefficients(cache,nset) ;
992
993 if (!cache._rangeProjList.empty()) {
994
995 for (std::size_t i = 0; i < _pdfList.size(); ++i) {
996 auto const& r1 = static_cast<RooAbsReal&>(cache._refRangeProjList[i]);
997 auto const& r2 = static_cast<RooAbsReal&>(cache._rangeProjList[i]);
998 double ncomp = _allExtendable ? static_cast<RooAbsPdf&>(_pdfList[i]).expectedEvents(nset)
999 : static_cast<RooAbsReal&>(_coefList[i]).getVal(nset);
1000 expectedTotal += (r2.getVal()/r1.getVal()) * ncomp ;
1001
1002 }
1003
1004 } else {
1005
1006 if (_allExtendable) {
1007 for(auto const& arg : _pdfList) {
1008 expectedTotal += static_cast<RooAbsPdf*>(arg)->expectedEvents(nset) ;
1009 }
1010 } else {
1011 for(auto const& arg : _coefList) {
1012 expectedTotal += static_cast<RooAbsReal*>(arg)->getVal(nset) ;
1013 }
1014 }
1015
1016 }
1017 return expectedTotal ;
1018}
1019
1020
1021
1022////////////////////////////////////////////////////////////////////////////////
1023/// Interface function used by test statistics to freeze choice of observables
1024/// for interpretation of fraction coefficients
1025
1026void RooAddPdf::selectNormalization(const RooArgSet* depSet, bool force)
1027{
1028
1029 if (!force && !_refCoefNorm.empty()) {
1030 return ;
1031 }
1032
1033 if (!depSet) {
1035 return ;
1036 }
1037
1038 fixCoefNormalization(*std::unique_ptr<RooArgSet>{getObservables(depSet)}) ;
1039}
1040
1041
1042
1043////////////////////////////////////////////////////////////////////////////////
1044/// Interface function used by test statistics to freeze choice of range
1045/// for interpretation of fraction coefficients
1046
1047void RooAddPdf::selectNormalizationRange(const char* rangeName, bool force)
1048{
1049 if (!force && _refCoefRangeName) {
1050 return ;
1051 }
1052
1053 fixCoefRange(rangeName) ;
1054}
1055
1056
1057
1058////////////////////////////////////////////////////////////////////////////////
1059/// Return specialized context to efficiently generate toy events from RooAddPdfs
1060/// return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ; // WVE DEBUG
1061
1063 const RooArgSet* auxProto, bool verbose) const
1064{
1065 return new RooAddGenContext(*this,vars,prototype,auxProto,verbose) ;
1066}
1067
1068
1069
1070////////////////////////////////////////////////////////////////////////////////
1071/// List all RooAbsArg derived contents in this cache element
1072
1074{
1075 RooArgList allNodes;
1076 allNodes.add(_projList) ;
1077 allNodes.add(_suppProjList) ;
1078 allNodes.add(_refRangeProjList) ;
1079 allNodes.add(_rangeProjList) ;
1080
1081 return allNodes ;
1082}
1083
1084
1085
1086////////////////////////////////////////////////////////////////////////////////
1087/// Loop over components for plot sampling hints and merge them if there are multiple
1088
1089std::list<Double_t>* RooAddPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
1090{
1091 std::unique_ptr<std::list<Double_t>> sumHint = nullptr ;
1092 bool needClean = false;
1093
1094 // Loop over components pdf
1095 for (const auto arg : _pdfList) {
1096 auto pdf = static_cast<const RooAbsPdf*>(arg);
1097
1098 std::unique_ptr<std::list<Double_t>> pdfHint{pdf->plotSamplingHint(obs,xlo,xhi)} ;
1099
1100 // Process hint
1101 if (pdfHint) {
1102 if (!sumHint) {
1103
1104 // If this is the first hint, then just save it
1105 sumHint = std::move(pdfHint) ;
1106
1107 } else {
1108
1109 auto newSumHint = std::make_unique<std::list<Double_t>>(sumHint->size()+pdfHint->size());
1110
1111 // Merge hints into temporary array
1112 merge(pdfHint->begin(),pdfHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
1113
1114 // Copy merged array without duplicates to new sumHintArrau
1115 sumHint = std::move(newSumHint) ;
1116 needClean = true ;
1117
1118 }
1119 }
1120 }
1121 if (needClean) {
1122 sumHint->erase(std::unique(sumHint->begin(),sumHint->end()), sumHint->end()) ;
1123 }
1124
1125 return sumHint.release() ;
1126}
1127
1128
1129////////////////////////////////////////////////////////////////////////////////
1130/// Loop over components for plot sampling hints and merge them if there are multiple
1131
1132std::list<Double_t>* RooAddPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
1133{
1134 std::unique_ptr<list<Double_t>> sumBinB = nullptr ;
1135 bool needClean = false;
1136
1137 // Loop over components pdf
1138 for (auto arg : _pdfList) {
1139 auto pdf = static_cast<const RooAbsPdf *>(arg);
1140 std::unique_ptr<list<Double_t>> pdfBinB{pdf->binBoundaries(obs,xlo,xhi)};
1141
1142 // Process hint
1143 if (pdfBinB) {
1144 if (!sumBinB) {
1145
1146 // If this is the first hint, then just save it
1147 sumBinB = std::move(pdfBinB) ;
1148
1149 } else {
1150
1151 auto newSumBinB = std::make_unique<list<Double_t>>(sumBinB->size()+pdfBinB->size()) ;
1152
1153 // Merge hints into temporary array
1154 merge(pdfBinB->begin(),pdfBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
1155
1156 // Copy merged array without duplicates to new sumBinBArrau
1157 sumBinB = std::move(newSumBinB) ;
1158 needClean = true ;
1159 }
1160 }
1161 }
1162
1163 // Remove consecutive duplicates
1164 if (needClean) {
1165 sumBinB->erase(std::unique(sumBinB->begin(),sumBinB->end()), sumBinB->end()) ;
1166 }
1167
1168 return sumBinB.release() ;
1169}
1170
1171
1172////////////////////////////////////////////////////////////////////////////////
1173/// If all components that depend on obs are binned, so is their sum.
1175{
1176 for (const auto arg : _pdfList) {
1177 auto pdf = static_cast<const RooAbsPdf*>(arg);
1178 if (pdf->dependsOn(obs) && !pdf->isBinnedDistribution(obs)) {
1179 return false ;
1180 }
1181 }
1182
1183 return true ;
1184}
1185
1186
1187////////////////////////////////////////////////////////////////////////////////
1188/// Label OK'ed components of a RooAddPdf with cache-and-track
1189
1191{
1192 RooFIter aiter = pdfList().fwdIterator() ;
1193 RooAbsArg* aarg ;
1194 while ((aarg=aiter.next())) {
1195 if (aarg->canNodeBeCached()==Always) {
1196 trackNodes.add(*aarg) ;
1197 //cout << "tracking node RooAddPdf component " << aarg->IsA()->GetName() << "::" << aarg->GetName() << endl ;
1198 }
1199 }
1200}
1201
1202
1203
1204////////////////////////////////////////////////////////////////////////////////
1205/// Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the
1206/// product operator construction
1207
1208void RooAddPdf::printMetaArgs(ostream& os) const
1209{
1210 bool first(true) ;
1211
1212 if (!_coefList.empty()) {
1213 for (std::size_t i = 0; i < _pdfList.size(); ++i ) {
1214 const RooAbsArg * coef = _coefList.at(i);
1215 const RooAbsArg * pdf = _pdfList.at(i);
1216 if (!first) {
1217 os << " + " ;
1218 } else {
1219 first = false ;
1220 }
1221
1222 if (i < _coefList.size()) {
1223 os << coef->GetName() << " * " << pdf->GetName();
1224 } else {
1225 os << "[%] * " << pdf->GetName();
1226 }
1227 }
1228 } else {
1229
1230 for (const auto pdf : _pdfList) {
1231 if (!first) {
1232 os << " + " ;
1233 } else {
1234 first = false ;
1235 }
1236 os << pdf->GetName() ;
1237 }
1238 }
1239
1240 os << " " ;
1241}
#define cxcoutD(a)
#define oocoutW(o, a)
#define coutW(a)
#define dologD(a)
#define coutE(a)
#define ccoutD(a)
#define TRACE_CREATE
Definition RooTrace.h:23
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
const std::vector< Int_t > & retrieve(Int_t masterCode) const
Retrieve the array of integer codes associated with the given master code.
Int_t store(const std::vector< Int_t > &codeList, RooArgSet *set1=0, RooArgSet *set2=0, RooArgSet *set3=0, RooArgSet *set4=0)
Store given arrays of integer codes, and up to four RooArgSets in the registry (each setX pointer may...
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
void clearValueAndShapeDirty() const
Definition RooAbsArg.h:611
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
virtual CacheMode canNodeBeCached() const
Definition RooAbsArg.h:427
friend class RooArgSet
Definition RooAbsArg.h:642
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
Bool_t isValueDirty() const
Definition RooAbsArg.h:436
void setOperMode(OperMode mode, Bool_t recurseADirty=kTRUE)
Set the operation mode of this node.
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:499
Storage_t const & get() const
Const access to the underlying stl container.
RooFIter fwdIterator() const
One-time forward iterator.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add an argument and transfer the ownership to the collection.
Storage_t::size_type size() const
RooAbsArg * first() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
const char * GetName() const
Returns name of object.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
TString _normRange
MC generator configuration specific for this object.
Definition RooAbsPdf.h:391
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:354
RooAbsReal * _norm
Definition RooAbsPdf.h:353
friend class RooRealIntegral
Definition RooAbsPdf.h:346
virtual Bool_t syncNormalization(const RooArgSet *dset, Bool_t adjustProxies=kTRUE) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
Int_t _errorCount
Definition RooAbsPdf.h:383
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
static Int_t _verboseEval
Definition RooAbsPdf.h:347
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
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 listed in ...
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
Double_t _value
Definition RooAbsReal.h:484
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
Transient cache with transformed values of coefficients.
Definition RooAddPdf.h:105
virtual RooArgList containedArgs(Action)
List all RooAbsArg derived contents in this cache element.
RooArgList _rangeProjList
Definition RooAddPdf.h:115
RooArgList _refRangeProjList
Definition RooAddPdf.h:114
RooArgList _projList
Definition RooAddPdf.h:112
RooArgList _suppNormList
Definition RooAddPdf.h:109
RooArgList _suppProjList
Definition RooAddPdf.h:113
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooListProxy _coefList
Definition RooAddPdf.h:141
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, bool verbose=false) const override
Return specialized context to efficiently generate toy events from RooAddPdfs return RooAbsPdf::genCo...
bool _allExtendable
Definition RooAddPdf.h:145
friend class RooAddGenContext
Definition RooAddPdf.h:125
RooAICRegistry _codeReg
Definition RooAddPdf.h:138
RooFit::UniqueId< RooArgSet >::Value_t _idOfLastUsedNormSet
!
Definition RooAddPdf.h:159
std::pair< const RooArgSet *, CacheElem * > getNormAndCache(const RooArgSet *nset) const
Look up projection cache and per-PDF norm sets.
std::unique_ptr< const RooArgSet > _copyOfLastNormSet
!
Definition RooAddPdf.h:160
Int_t _coefErrCount
Definition RooAddPdf.h:148
bool _haveLastCoef
List of supplemental normalization factors.
Definition RooAddPdf.h:144
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=0) const override
Determine which part (if any) of given integral can be performed analytically.
void updateCoefficients(CacheElem &cache, const RooArgSet *nset) const
Update the coefficient values in the given cache element: calculate new remainder fraction,...
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the produ...
CacheElem * getProjCache(const RooArgSet *nset, const RooArgSet *iset=0, const char *rangeName=0) const
Manager of cache with coefficient projections and transformations.
void selectNormalization(const RooArgSet *depSet=0, bool force=false) override
Interface function used by test statistics to freeze choice of observables for interpretation of frac...
void finalizeConstruction()
void setCacheAndTrackHints(RooArgSet &) override
Label OK'ed components of a RooAddPdf with cache-and-track.
bool _recursive
Definition RooAddPdf.h:146
RooObjCacheManager _projCacheMgr
Definition RooAddPdf.h:120
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const override
Return analytical integral defined by given scenario code.
std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const override
Loop over components for plot sampling hints and merge them if there are multiple.
bool checkObservables(const RooArgSet *nset) const override
Check if PDF is valid for given normalization set.
void selectNormalizationRange(const char *rangeName=0, bool force=false) override
Interface function used by test statistics to freeze choice of range for interpretation of fraction c...
void fixCoefNormalization(const RooArgSet &refCoefNorm)
By default the interpretation of the fraction coefficients is performed in the contextual choice of o...
RooSetProxy _refCoefNorm
Definition RooAddPdf.h:98
void resetErrorCounters(Int_t resetValue=10) override
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
double getValV(const RooArgSet *set=nullptr) const override
Calculate and return the current value.
void fixCoefRange(const char *rangeName)
By default, fraction coefficients are assumed to refer to the default fit range.
RooListProxy _pdfList
Registry of component analytical integration codes.
Definition RooAddPdf.h:140
TNamed * _refCoefRangeName
Definition RooAddPdf.h:99
const RooArgList & pdfList() const
Definition RooAddPdf.h:67
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute addition of PDFs in batches.
bool isBinnedDistribution(const RooArgSet &obs) const override
If all components that depend on obs are binned, so is their sum.
std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const override
Loop over components for plot sampling hints and merge them if there are multiple.
Double_t expectedEvents(const RooArgSet *nset) const override
Return expected number of events for extended likelihood calculation, which is the sum of all coeffic...
bool _projectCoefs
Definition RooAddPdf.h:101
std::vector< double > _coefCache
Definition RooAddPdf.h:102
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
RooFit::UniqueId< RooArgSet > const & uniqueId() const
Returns a unique ID that is different for every instantiated RooArgSet.
Definition RooArgSet.h:171
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
void reset()
Clear the cache.
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:36
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
static RooMsgService & instance()
Return reference to singleton instance.
static Int_t _debugCount
Bool_t isActive(const RooAbsArg *self, RooFit::MsgTopic facility, RooFit::MsgLevel level)
Check if logging is active for given object/topic/RooFitMsgLevel combination.
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
static RooConstVar & value(Double_t value)
Return a constant value object with given value.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
Class RooRecursiveFraction is a RooAbsReal implementation that calculates the plain fraction of sum o...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Ssiz_t Length() const
Definition TString.h:410
const char * Data() const
Definition TString.h:369
RooConstVar & RooConst(Double_t val)
const Int_t n
Definition legend1.C:16
std::vector< RooSpan< const double > > VarVector
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
std::vector< double > ArgVector
Definition first.py:1
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output(int code)
Definition gifencode.c:226