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