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