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