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