Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
22 RooAddPdf 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 
28 or
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 
33 The first form is for extended likelihood fits, where the
34 expected number of events is \f$ \sum_i c_i \f$. The coefficients \f$ c_i \f$
35 can either be explicitly provided, or, if all components support
36 extended likelihood fits, they can be calculated the contribution
37 of each PDF to the total number of expected events.
38 
39 In the second form, the sum of the coefficients is required to be one or less,
40 and the coefficient of the last PDF is calculated from that condition.
41 
42 It is also possible to parameterize the coefficients recursively
43 
44 \f[
45 c_1*\mathrm{PDF}_1 + (1-c_1)(c_2*\mathrm{PDF}_2 + (1-c_2)*(c_3*\mathrm{PDF}_3 + ...))
46 \f]
47 
48 In this form the sum of the coefficients is always less than 1.0
49 for all possible values of the individual coefficients between 0 and 1.
50 
51 RooAddPdf relies on each component PDF to be normalized and will perform
52 no normalization other than calculating the proper last coefficient \f$ c_n \f$, if requested.
53 An (enforced) condition for this assumption is that each \f$ \mathrm{PDF}_i \f$ is independent
54 of 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"
80 #include "RooRecursiveFraction.h"
81 #include "RooGlobalFunc.h"
82 #include "RooRealIntegral.h"
83 #include "RooTrace.h"
84 
85 #include "Riostream.h"
86 #include <algorithm>
87 
88 
89 using namespace std;
90 
91 ClassImp(RooAddPdf);
92 ;
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Default constructor used for persistence
97 
98 RooAddPdf::RooAddPdf() :
99  _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
100  _refCoefRangeName(0),
101  _codeReg(10),
102  _snormList(0),
103  _recursive(kFALSE)
104 {
105  _coefErrCount = _errorCount ;
106  TRACE_CREATE
107 }
108 
109 
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 /// Dummy constructor
113 
114 RooAddPdf::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 {
128  _coefErrCount = _errorCount ;
129  TRACE_CREATE
130 }
131 
132 
133 
134 ////////////////////////////////////////////////////////////////////////////////
135 /// Constructor with two PDFs and one coefficient
136 
137 RooAddPdf::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());
156  _coefErrCount = _errorCount ;
157  TRACE_CREATE
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 
173 RooAddPdf::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 
200  Bool_t first(kTRUE) ;
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
262  _haveLastCoef=kTRUE ;
263  }
264 
265  } else {
266  _haveLastCoef=kTRUE ;
267  }
268 
269 
270  _coefCache.resize(_pdfList.size());
271  _coefErrCount = _errorCount ;
272  _recursive = recursiveFractions ;
273 
274  TRACE_CREATE
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 
285 RooAddPdf::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());
313  _coefErrCount = _errorCount ;
314  TRACE_CREATE
315 }
316 
317 
318 
319 
320 ////////////////////////////////////////////////////////////////////////////////
321 /// Copy constructor
322 
323 RooAddPdf::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());
337  _coefErrCount = _errorCount ;
338  TRACE_CREATE
339 }
340 
341 
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// Destructor
345 
346 RooAddPdf::~RooAddPdf()
347 {
348  TRACE_DESTROY
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 
363 void RooAddPdf::fixCoefNormalization(const RooArgSet& refCoefNorm)
364 {
365  if (refCoefNorm.getSize()==0) {
366  _projectCoefs = kFALSE ;
367  return ;
368  }
369  _projectCoefs = kTRUE ;
370 
371  _refCoefNorm.removeAll() ;
372  _refCoefNorm.add(refCoefNorm) ;
373 
374  _projCacheMgr.reset() ;
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 
390 void RooAddPdf::fixCoefRange(const char* rangeName)
391 {
392  _refCoefRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
393  if (_refCoefRangeName) _projectCoefs = kTRUE ;
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 
409 RooAddPdf::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) {
502  rangeName = RooNameReg::str(_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 
652 void 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  {
744  RooAbsReal::GlobalSelectComponentRAII compRAII(true);
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 
768  if ((RooMsgService::_debugCount>0) && RooMsgService::instance().isActive(this,RooFit::Caching,RooFit::DEBUG)) {
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 
787 std::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 
807 Double_t RooAddPdf::evaluate() const
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 
839 RooSpan<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 
873 void RooAddPdf::resetErrorCounters(Int_t resetValue)
874 {
875  RooAbsPdf::resetErrorCounters(resetValue) ;
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 
886 Bool_t RooAddPdf::checkObservables(const RooArgSet* nset) const
887 {
888  Bool_t ret(kFALSE) ;
889 
890  for (int i = 0; i < _pdfList.getSize(); ++i) {
891  auto pdf = static_cast<const RooAbsPdf *>(_pdfList.at(i));
892  auto coef = static_cast<const RooAbsReal*>(_coefList.at(i));
893  if (pdf->observableOverlaps(nset,*coef)) {
894  coutE(InputArguments) << "RooAddPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName()
895  << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
896  ret = kTRUE ;
897  }
898  }
899 
900  return ret ;
901 }
902 
903 
904 ////////////////////////////////////////////////////////////////////////////////
905 /// Determine which part (if any) of given integral can be performed analytically.
906 /// If any analytical integration is possible, return integration scenario code
907 ///
908 /// RooAddPdf queries each component PDF for its analytical integration capability of the requested
909 /// set ('allVars'). It finds the largest common set of variables that can be integrated
910 /// by all components. If such a set exists, it reconfirms that each component is capable of
911 /// analytically integrating the common set, and combines the components individual integration
912 /// codes into a single integration code valid for RooAddPdf.
913 
914 Int_t RooAddPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars,
915  const RooArgSet* normSet, const char* rangeName) const
916 {
917 
918  RooArgSet* allDepVars = getObservables(allVars) ;
919  RooArgSet allAnalVars(*allDepVars) ;
920  delete allDepVars ;
921 
922  Int_t n(0) ;
923 
924  // First iteration, determine what each component can integrate analytically
925  for (const auto pdfArg : _pdfList) {
926  auto pdf = static_cast<const RooAbsPdf *>(pdfArg);
927  RooArgSet subAnalVars ;
928  pdf->getAnalyticalIntegralWN(allVars,subAnalVars,normSet,rangeName) ;
929 
930  // Observables that cannot be integrated analytically by this component are dropped from the common list
931  for (const auto arg : allVars) {
932  if (!subAnalVars.find(arg->GetName()) && pdf->dependsOn(*arg)) {
933  allAnalVars.remove(*arg,kTRUE,kTRUE) ;
934  }
935  }
936  n++ ;
937  }
938 
939  // If no observables can be integrated analytically, return code 0 here
940  if (allAnalVars.getSize()==0) {
941  return 0 ;
942  }
943 
944 
945  // Now retrieve codes for integration over common set of analytically integrable observables for each component
946  n=0 ;
947  std::vector<Int_t> subCode(_pdfList.getSize());
948  Bool_t allOK(kTRUE) ;
949  for (const auto arg : _pdfList) {
950  auto pdf = static_cast<const RooAbsPdf *>(arg);
951  RooArgSet subAnalVars ;
952  RooArgSet* allAnalVars2 = pdf->getObservables(allAnalVars) ;
953  subCode[n] = pdf->getAnalyticalIntegralWN(*allAnalVars2,subAnalVars,normSet,rangeName) ;
954  if (subCode[n]==0 && allAnalVars2->getSize()>0) {
955  coutE(InputArguments) << "RooAddPdf::getAnalyticalIntegral(" << GetName() << ") WARNING: component PDF " << pdf->GetName()
956  << " advertises inconsistent set of integrals (e.g. (X,Y) but not X or Y individually."
957  << " Distributed analytical integration disabled. Please fix PDF" << endl ;
958  allOK = kFALSE ;
959  }
960  delete allAnalVars2 ;
961  n++ ;
962  }
963  if (!allOK) {
964  return 0 ;
965  }
966 
967  // Mare all analytically integrated observables as such
968  analVars.add(allAnalVars) ;
969 
970  // Store set of variables analytically integrated
971  RooArgSet* intSet = new RooArgSet(allAnalVars) ;
972  Int_t masterCode = _codeReg.store(subCode,intSet)+1 ;
973 
974  return masterCode ;
975 }
976 
977 
978 
979 ////////////////////////////////////////////////////////////////////////////////
980 /// Return analytical integral defined by given scenario code
981 
982 Double_t RooAddPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
983 {
984  // WVE needs adaptation to handle new rangeName feature
985  if (code==0) {
986  return getVal(normSet) ;
987  }
988 
989  // Retrieve analytical integration subCodes and set of observabels integrated over
990  RooArgSet* intSet ;
991  const std::vector<Int_t>& subCode = _codeReg.retrieve(code-1,intSet) ;
992  if (subCode.empty()) {
993  coutE(InputArguments) << "RooAddPdf::analyticalIntegral(" << GetName() << "): ERROR unrecognized integration code, " << code << endl ;
994  assert(0) ;
995  }
996 
997  cxcoutD(Caching) << "RooAddPdf::aiWN(" << GetName() << ") calling getProjCache with nset = " << (normSet?*normSet:RooArgSet()) << endl ;
998 
999  if ((normSet==0 || normSet->getSize()==0) && _refCoefNorm.getSize()>0) {
1000 // cout << "WVE integration of RooAddPdf without normalization, but have reference set, using ref set for normalization" << endl ;
1001  normSet = &_refCoefNorm ;
1002  }
1003 
1004  CacheElem* cache = getProjCache(normSet,intSet,0) ; // WVE rangename here?
1005  updateCoefficients(*cache,normSet) ;
1006 
1007  // Calculate the current value of this object
1008  Double_t value(0) ;
1009 
1010  // Do running sum of coef/pdf pairs, calculate lastCoef.
1011  Double_t snormVal ;
1012 
1013  //cout << "ROP::aIWN updateCoefCache with rangeName = " << (rangeName?rangeName:"<null>") << endl ;
1014  RooArgList* snormSet = (cache->_suppNormList.getSize()>0) ? &cache->_suppNormList : 0 ;
1015  for (int i = 0; i < _pdfList.getSize(); ++i ) {
1016  auto pdf = static_cast<const RooAbsPdf*>(_pdfList.at(i));
1017 
1018  if (_coefCache[i]) {
1019  snormVal = snormSet ? ((RooAbsReal*) cache->_suppNormList.at(i))->getVal() : 1.0 ;
1020 
1021  // WVE swap this?
1022  Double_t val = pdf->analyticalIntegralWN(subCode[i],normSet,rangeName) ;
1023  if (pdf->isSelectedComp()) {
1024  value += val*_coefCache[i]/snormVal ;
1025  }
1026  }
1027  }
1028 
1029  return value ;
1030 }
1031 
1032 
1033 
1034 ////////////////////////////////////////////////////////////////////////////////
1035 /// Return the number of expected events, which is either the sum of all coefficients
1036 /// or the sum of the components extended terms, multiplied with the fraction that
1037 /// is in the current range w.r.t the reference range
1038 
1039 Double_t RooAddPdf::expectedEvents(const RooArgSet* nset) const
1040 {
1041  Double_t expectedTotal(0.0);
1042 
1043  cxcoutD(Caching) << "RooAddPdf::expectedEvents(" << GetName() << ") calling getProjCache with nset = " << (nset?*nset:RooArgSet()) << endl ;
1044  CacheElem* cache = getProjCache(nset) ;
1045  updateCoefficients(*cache,nset) ;
1046 
1047  if (cache->_rangeProjList.getSize()>0) {
1048 
1049  RooFIter iter1 = cache->_refRangeProjList.fwdIterator() ;
1050  RooFIter iter2 = cache->_rangeProjList.fwdIterator() ;
1051  RooFIter iter3 = _pdfList.fwdIterator() ;
1052 
1053  if (_allExtendable) {
1054 
1055  RooAbsPdf* pdf ;
1056  while ((pdf=(RooAbsPdf*)iter3.next())) {
1057  RooAbsReal* r1 = (RooAbsReal*)iter1.next() ;
1058  RooAbsReal* r2 = (RooAbsReal*)iter2.next() ;
1059  expectedTotal += (r2->getVal()/r1->getVal()) * pdf->expectedEvents(nset) ;
1060  }
1061 
1062  } else {
1063 
1064  RooFIter citer = _coefList.fwdIterator() ;
1065  RooAbsReal* coef ;
1066  while((coef=(RooAbsReal*)citer.next())) {
1067  Double_t ncomp = coef->getVal(nset) ;
1068  RooAbsReal* r1 = (RooAbsReal*)iter1.next() ;
1069  RooAbsReal* r2 = (RooAbsReal*)iter2.next() ;
1070  expectedTotal += (r2->getVal()/r1->getVal()) * ncomp ;
1071  }
1072 
1073  }
1074 
1075 
1076 
1077  } else {
1078 
1079  if (_allExtendable) {
1080 
1081  RooFIter iter = _pdfList.fwdIterator() ;
1082  RooAbsPdf* pdf ;
1083  while((pdf=(RooAbsPdf*)iter.next())) {
1084  expectedTotal += pdf->expectedEvents(nset) ;
1085  }
1086 
1087  } else {
1088 
1089  RooFIter citer = _coefList.fwdIterator() ;
1090  RooAbsReal* coef ;
1091  while((coef=(RooAbsReal*)citer.next())) {
1092  Double_t ncomp = coef->getVal(nset) ;
1093  expectedTotal += ncomp ;
1094  }
1095 
1096  }
1097 
1098  }
1099  return expectedTotal ;
1100 }
1101 
1102 
1103 
1104 ////////////////////////////////////////////////////////////////////////////////
1105 /// Interface function used by test statistics to freeze choice of observables
1106 /// for interpretation of fraction coefficients
1107 
1108 void RooAddPdf::selectNormalization(const RooArgSet* depSet, Bool_t force)
1109 {
1110 
1111  if (!force && _refCoefNorm.getSize()!=0) {
1112  return ;
1113  }
1114 
1115  if (!depSet) {
1116  fixCoefNormalization(RooArgSet()) ;
1117  return ;
1118  }
1119 
1120  RooArgSet* myDepSet = getObservables(depSet) ;
1121  fixCoefNormalization(*myDepSet) ;
1122  delete myDepSet ;
1123 }
1124 
1125 
1126 
1127 ////////////////////////////////////////////////////////////////////////////////
1128 /// Interface function used by test statistics to freeze choice of range
1129 /// for interpretation of fraction coefficients
1130 
1131 void RooAddPdf::selectNormalizationRange(const char* rangeName, Bool_t force)
1132 {
1133  if (!force && _refCoefRangeName) {
1134  return ;
1135  }
1136 
1137  fixCoefRange(rangeName) ;
1138 }
1139 
1140 
1141 
1142 ////////////////////////////////////////////////////////////////////////////////
1143 /// Return specialized context to efficiently generate toy events from RooAddPdfs
1144 /// return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ; // WVE DEBUG
1145 
1146 RooAbsGenContext* RooAddPdf::genContext(const RooArgSet &vars, const RooDataSet *prototype,
1147  const RooArgSet* auxProto, Bool_t verbose) const
1148 {
1149  return new RooAddGenContext(*this,vars,prototype,auxProto,verbose) ;
1150 }
1151 
1152 
1153 
1154 ////////////////////////////////////////////////////////////////////////////////
1155 /// List all RooAbsArg derived contents in this cache element
1156 
1157 RooArgList RooAddPdf::CacheElem::containedArgs(Action)
1158 {
1159  RooArgList allNodes;
1160  allNodes.add(_projList) ;
1161  allNodes.add(_suppProjList) ;
1162  allNodes.add(_refRangeProjList) ;
1163  allNodes.add(_rangeProjList) ;
1164 
1165  return allNodes ;
1166 }
1167 
1168 
1169 
1170 ////////////////////////////////////////////////////////////////////////////////
1171 /// Loop over components for plot sampling hints and merge them if there are multiple
1172 
1173 std::list<Double_t>* RooAddPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
1174 {
1175  list<Double_t>* sumHint = 0 ;
1176  Bool_t needClean(kFALSE) ;
1177 
1178  // Loop over components pdf
1179  for (const auto arg : _pdfList) {
1180  auto pdf = static_cast<const RooAbsPdf*>(arg);
1181 
1182  list<Double_t>* pdfHint = pdf->plotSamplingHint(obs,xlo,xhi) ;
1183 
1184  // Process hint
1185  if (pdfHint) {
1186  if (!sumHint) {
1187 
1188  // If this is the first hint, then just save it
1189  sumHint = pdfHint ;
1190 
1191  } else {
1192 
1193  list<Double_t>* newSumHint = new list<Double_t>(sumHint->size()+pdfHint->size()) ;
1194 
1195  // Merge hints into temporary array
1196  merge(pdfHint->begin(),pdfHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
1197 
1198  // Copy merged array without duplicates to new sumHintArrau
1199  delete sumHint ;
1200  sumHint = newSumHint ;
1201  needClean = kTRUE ;
1202 
1203  }
1204  }
1205  }
1206  if (needClean) {
1207  list<Double_t>::iterator new_end = unique(sumHint->begin(),sumHint->end()) ;
1208  sumHint->erase(new_end,sumHint->end()) ;
1209  }
1210 
1211  return sumHint ;
1212 }
1213 
1214 
1215 ////////////////////////////////////////////////////////////////////////////////
1216 /// Loop over components for plot sampling hints and merge them if there are multiple
1217 
1218 std::list<Double_t>* RooAddPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
1219 {
1220  list<Double_t>* sumBinB = 0 ;
1221  Bool_t needClean(kFALSE) ;
1222 
1223  // Loop over components pdf
1224  for (auto arg : _pdfList) {
1225  auto pdf = static_cast<const RooAbsPdf *>(arg);
1226  list<Double_t>* pdfBinB = pdf->binBoundaries(obs,xlo,xhi) ;
1227 
1228  // Process hint
1229  if (pdfBinB) {
1230  if (!sumBinB) {
1231 
1232  // If this is the first hint, then just save it
1233  sumBinB = pdfBinB ;
1234 
1235  } else {
1236 
1237  list<Double_t>* newSumBinB = new list<Double_t>(sumBinB->size()+pdfBinB->size()) ;
1238 
1239  // Merge hints into temporary array
1240  merge(pdfBinB->begin(),pdfBinB->end(),sumBinB->begin(),sumBinB->end(),newSumBinB->begin()) ;
1241 
1242  // Copy merged array without duplicates to new sumBinBArrau
1243  delete sumBinB ;
1244  delete pdfBinB ;
1245  sumBinB = newSumBinB ;
1246  needClean = kTRUE ;
1247  }
1248  }
1249  }
1250 
1251  // Remove consecutive duplicates
1252  if (needClean) {
1253  list<Double_t>::iterator new_end = unique(sumBinB->begin(),sumBinB->end()) ;
1254  sumBinB->erase(new_end,sumBinB->end()) ;
1255  }
1256 
1257  return sumBinB ;
1258 }
1259 
1260 
1261 ////////////////////////////////////////////////////////////////////////////////
1262 /// If all components that depend on obs are binned that so is the product
1263 
1264 Bool_t RooAddPdf::isBinnedDistribution(const RooArgSet& obs) const
1265 {
1266  for (const auto arg : _pdfList) {
1267  auto pdf = static_cast<const RooAbsPdf*>(arg);
1268  if (pdf->dependsOn(obs) && !pdf->isBinnedDistribution(obs)) {
1269  return kFALSE ;
1270  }
1271  }
1272 
1273  return kTRUE ;
1274 }
1275 
1276 
1277 ////////////////////////////////////////////////////////////////////////////////
1278 /// Label OK'ed components of a RooAddPdf with cache-and-track
1279 
1280 void RooAddPdf::setCacheAndTrackHints(RooArgSet& trackNodes)
1281 {
1282  RooFIter aiter = pdfList().fwdIterator() ;
1283  RooAbsArg* aarg ;
1284  while ((aarg=aiter.next())) {
1285  if (aarg->canNodeBeCached()==Always) {
1286  trackNodes.add(*aarg) ;
1287  //cout << "tracking node RooAddPdf component " << aarg->IsA()->GetName() << "::" << aarg->GetName() << endl ;
1288  }
1289  }
1290 }
1291 
1292 
1293 
1294 ////////////////////////////////////////////////////////////////////////////////
1295 /// Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the
1296 /// product operator construction
1297 
1298 void RooAddPdf::printMetaArgs(ostream& os) const
1299 {
1300  Bool_t first(kTRUE) ;
1301 
1302  if (_coefList.getSize() != 0) {
1303  for (int i = 0; i < _pdfList.getSize(); ++i ) {
1304  const RooAbsArg * coef = _coefList.at(i);
1305  const RooAbsArg * pdf = _pdfList.at(i);
1306  if (!first) {
1307  os << " + " ;
1308  } else {
1309  first = kFALSE ;
1310  }
1311 
1312  if (i < _coefList.getSize()) {
1313  os << coef->GetName() << " * " << pdf->GetName();
1314  } else {
1315  os << "[%] * " << pdf->GetName();
1316  }
1317  }
1318  } else {
1319 
1320  for (const auto pdf : _pdfList) {
1321  if (!first) {
1322  os << " + " ;
1323  } else {
1324  first = kFALSE ;
1325  }
1326  os << pdf->GetName() ;
1327  }
1328  }
1329 
1330  os << " " ;
1331 }