Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooAddModel.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 // RooAddModel 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 // RooAddPdf relies on each component PDF to be normalized and will perform
37 // no normalization other than calculating the proper last coefficient c_n, if requested.
38 // An (enforced) condition for this assuption is that each PDF_i is independent
39 // of each coefficient_i.
40 //
41 //
42 
43 #include "RooFit.h"
44 #include "RooMsgService.h"
45 
46 #include "TIterator.h"
47 #include "TIterator.h"
48 #include "TList.h"
49 #include "RooAddModel.h"
50 #include "RooDataSet.h"
51 #include "RooRealProxy.h"
52 #include "RooPlot.h"
53 #include "RooRealVar.h"
54 #include "RooAddGenContext.h"
55 #include "RooRealConstant.h"
56 #include "RooNameReg.h"
57 #include "RooMsgService.h"
58 #include "RooRealIntegral.h"
59 
60 #include "Riostream.h"
61 
62 
63 
64 using namespace std;
65 
66 ClassImp(RooAddModel);
67 ;
68 
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 
72 RooAddModel::RooAddModel() :
73  _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
74  _refCoefRangeName(0),
75  _codeReg(10),
76  _snormList(0)
77 {
78  _pdfIter = _pdfList.createIterator() ;
79  _coefIter = _coefList.createIterator() ;
80 
81  _coefCache = new Double_t[10] ;
82  _coefErrCount = _errorCount ;
83 }
84 
85 
86 
87 ////////////////////////////////////////////////////////////////////////////////
88 /// Generic constructor from list of PDFs and list of coefficients.
89 /// Each pdf list element (i) is paired with coefficient list element (i).
90 /// The number of coefficients must be either equal to the number of PDFs,
91 /// in which case extended MLL fitting is enabled, or be one less.
92 ///
93 /// All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal
94 
95 RooAddModel::RooAddModel(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, Bool_t ownPdfList) :
96  RooResolutionModel(name,title,((RooResolutionModel*)inPdfList.at(0))->convVar()),
97  _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
98  _refCoefRangeName(0),
99  _projectCoefs(kFALSE),
100  _projCacheMgr(this,10),
101  _intCacheMgr(this,10),
102  _codeReg(10),
103  _pdfList("!pdfs","List of PDFs",this),
104  _coefList("!coefficients","List of coefficients",this),
105  _haveLastCoef(kFALSE),
106  _allExtendable(kFALSE)
107 {
108  if (inPdfList.getSize()>inCoefList.getSize()+1) {
109  coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName()
110  << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
111  assert(0) ;
112  }
113 
114  _pdfIter = _pdfList.createIterator() ;
115  _coefIter = _coefList.createIterator() ;
116 
117  // Constructor with N PDFs and N or N-1 coefs
118  TIterator* pdfIter = inPdfList.createIterator() ;
119  TIterator* coefIter = inCoefList.createIterator() ;
120  RooAbsPdf* pdf ;
121  RooAbsReal* coef ;
122 
123  while((coef = (RooAbsPdf*)coefIter->Next())) {
124  pdf = (RooAbsPdf*) pdfIter->Next() ;
125  if (!pdf) {
126  coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName()
127  << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
128  assert(0) ;
129  }
130  if (!dynamic_cast<RooAbsReal*>(coef)) {
131  coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal, ignored" << endl ;
132  continue ;
133  }
134  if (!dynamic_cast<RooAbsReal*>(pdf)) {
135  coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") pdf " << pdf->GetName() << " is not of type RooAbsPdf, ignored" << endl ;
136  continue ;
137  }
138  _pdfList.add(*pdf) ;
139  _coefList.add(*coef) ;
140  }
141 
142  pdf = (RooAbsPdf*) pdfIter->Next() ;
143  if (pdf) {
144  if (!dynamic_cast<RooAbsReal*>(pdf)) {
145  coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") last pdf " << coef->GetName() << " is not of type RooAbsPdf, fatal error" << endl ;
146  assert(0) ;
147  }
148  _pdfList.add(*pdf) ;
149  } else {
150  _haveLastCoef=kTRUE ;
151  }
152 
153  delete pdfIter ;
154  delete coefIter ;
155 
156  _coefCache = new Double_t[_pdfList.getSize()] ;
157  _coefErrCount = _errorCount ;
158 
159  if (ownPdfList) {
160  _ownedComps.addOwned(_pdfList) ;
161  }
162 
163 }
164 
165 
166 
167 ////////////////////////////////////////////////////////////////////////////////
168 /// Copy constructor
169 
170 RooAddModel::RooAddModel(const RooAddModel& other, const char* name) :
171  RooResolutionModel(other,name),
172  _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
173  _refCoefRangeName((TNamed*)other._refCoefRangeName),
174  _projectCoefs(other._projectCoefs),
175  _projCacheMgr(other._projCacheMgr,this),
176  _intCacheMgr(other._intCacheMgr,this),
177  _codeReg(other._codeReg),
178  _pdfList("!pdfs",this,other._pdfList),
179  _coefList("!coefficients",this,other._coefList),
180  _haveLastCoef(other._haveLastCoef),
181  _allExtendable(other._allExtendable)
182 {
183  _pdfIter = _pdfList.createIterator() ;
184  _coefIter = _coefList.createIterator() ;
185  _coefCache = new Double_t[_pdfList.getSize()] ;
186  _coefErrCount = _errorCount ;
187 }
188 
189 
190 
191 ////////////////////////////////////////////////////////////////////////////////
192 /// Destructor
193 
194 RooAddModel::~RooAddModel()
195 {
196  delete _pdfIter ;
197  delete _coefIter ;
198 
199  if (_coefCache) delete[] _coefCache ;
200 }
201 
202 
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 /// By default the interpretation of the fraction coefficients is
206 /// performed in the contextual choice of observables. This makes the
207 /// shape of the p.d.f explicitly dependent on the choice of
208 /// observables. This method instructs RooAddPdf to freeze the
209 /// interpretation of the coefficients to be done in the given set of
210 /// observables. If frozen, fractions are automatically transformed
211 /// from the reference normalization set to the contextual normalization
212 /// set by ratios of integrals
213 
214 void RooAddModel::fixCoefNormalization(const RooArgSet& refCoefNorm)
215 {
216  if (refCoefNorm.getSize()==0) {
217  _projectCoefs = kFALSE ;
218  return ;
219  }
220  _projectCoefs = kTRUE ;
221 
222  _refCoefNorm.removeAll() ;
223  _refCoefNorm.add(refCoefNorm) ;
224 
225  _projCacheMgr.reset() ;
226 }
227 
228 
229 
230 ////////////////////////////////////////////////////////////////////////////////
231 /// By default the interpretation of the fraction coefficients is
232 /// performed in the default range. This make the shape of a RooAddPdf
233 /// explicitly dependent on the range of the observables. To allow
234 /// a range independent definition of the fraction this function
235 /// instructs RooAddPdf to freeze its interpretation in the given
236 /// named range. If the current normalization range is different
237 /// from the reference range, the appropriate fraction coefficients
238 /// are automically calculation from the reference fractions using
239 /// ratios if integrals
240 
241 void RooAddModel::fixCoefRange(const char* rangeName)
242 {
243  _refCoefRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
244  if (_refCoefRangeName) _projectCoefs = kTRUE ;
245 }
246 
247 
248 
249 ////////////////////////////////////////////////////////////////////////////////
250 /// Instantiate a clone of this resolution model representing a convolution with given
251 /// basis function. The owners object name is incorporated in the clones name
252 /// to avoid multiple convolution objects with the same name in complex PDF structures.
253 ///
254 /// RooAddModel will clone all the component models to create a composite convolution object
255 
256 RooResolutionModel* RooAddModel::convolution(RooFormulaVar* inBasis, RooAbsArg* owner) const
257 {
258  // Check that primary variable of basis functions is our convolution variable
259  if (inBasis->getParameter(0) != x.absArg()) {
260  coutE(InputArguments) << "RooAddModel::convolution(" << GetName()
261  << ") convolution parameter of basis function and PDF don't match" << endl ;
262  ccoutE(InputArguments) << "basis->findServer(0) = " << inBasis->findServer(0) << " " << inBasis->findServer(0)->GetName() << endl ;
263  ccoutE(InputArguments) << "x.absArg() = " << x.absArg() << " " << x.absArg()->GetName() << endl ;
264  inBasis->Print("v") ;
265  return 0 ;
266  }
267 
268  TString newName(GetName()) ;
269  newName.Append("_conv_") ;
270  newName.Append(inBasis->GetName()) ;
271  newName.Append("_[") ;
272  newName.Append(owner->GetName()) ;
273  newName.Append("]") ;
274 
275  TString newTitle(GetTitle()) ;
276  newTitle.Append(" convoluted with basis function ") ;
277  newTitle.Append(inBasis->GetName()) ;
278 
279  _pdfIter->Reset() ;
280  RooResolutionModel* model ;
281  RooArgList modelList ;
282  while((model = (RooResolutionModel*)_pdfIter->Next())) {
283  // Create component convolution
284  RooResolutionModel* conv = model->convolution(inBasis,owner) ;
285  modelList.add(*conv) ;
286  }
287 
288  _coefIter->Reset() ;
289  RooAbsReal* coef ;
290  RooArgList theCoefList ;
291  while((coef = (RooAbsReal*)_coefIter->Next())) {
292  theCoefList.add(*coef) ;
293  }
294 
295  RooAddModel* convSum = new RooAddModel(newName,newTitle,modelList,theCoefList,kTRUE) ;
296  for (std::set<std::string>::const_iterator attrIt = _boolAttrib.begin();
297  attrIt != _boolAttrib.end(); ++attrIt) {
298  convSum->setAttribute((*attrIt).c_str()) ;
299  }
300  for (std::map<std::string,std::string>::const_iterator attrIt = _stringAttrib.begin();
301  attrIt != _stringAttrib.end(); ++attrIt) {
302  convSum->setStringAttribute((attrIt->first).c_str(), (attrIt->second).c_str()) ;
303  }
304  convSum->changeBasis(inBasis) ;
305  return convSum ;
306 }
307 
308 
309 
310 ////////////////////////////////////////////////////////////////////////////////
311 /// Return code for basis function representing by 'name' string.
312 /// The basis code of the first component model will be returned,
313 /// if the basis is supported by all components. Otherwise 0
314 /// is returned
315 
316 Int_t RooAddModel::basisCode(const char* name) const
317 {
318  TIterator* mIter = _pdfList.createIterator() ;
319  RooResolutionModel* model ;
320  Bool_t first(kTRUE), code(0) ;
321  while((model = (RooResolutionModel*)mIter->Next())) {
322  Int_t subCode = model->basisCode(name) ;
323  if (first) {
324  code = subCode ;
325  first = kFALSE ;
326  } else if (subCode==0) {
327  code = 0 ;
328  }
329  }
330  delete mIter ;
331 
332  return code ;
333 }
334 
335 
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated over iset
339 /// in range 'rangeName'. If cache element does not exist, create and fill it on the fly. The cache contains
340 /// suplemental normalization terms (in case not all added p.d.f.s have the same observables), projection
341 /// integrals to calculated transformed fraction coefficients when a frozen reference frame is provided
342 /// and projection integrals for similar transformations when a frozen reference range is provided.
343 
344 RooAddModel::CacheElem* RooAddModel::getProjCache(const RooArgSet* nset, const RooArgSet* iset, const char* rangeName) const
345 {
346  // Check if cache already exists
347  CacheElem* cache = (CacheElem*) _projCacheMgr.getObj(nset,iset,0,RooNameReg::ptr(rangeName)) ;
348  if (cache) {
349  return cache ;
350  }
351 
352  //Create new cache
353  cache = new CacheElem ;
354 
355  // *** PART 1 : Create supplemental normalization list ***
356 
357  // Retrieve the combined set of dependents of this PDF ;
358  RooArgSet *fullDepList = getObservables(nset) ;
359  if (iset) {
360  fullDepList->remove(*iset,kTRUE,kTRUE) ;
361  }
362 
363  // Fill with dummy unit RRVs for now
364  _pdfIter->Reset() ;
365  _coefIter->Reset() ;
366  RooAbsPdf* pdf ;
367  RooAbsReal* coef ;
368  while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
369  coef=(RooAbsPdf*)_coefIter->Next() ;
370 
371  // Start with full list of dependents
372  RooArgSet supNSet(*fullDepList) ;
373 
374  // Remove PDF dependents
375  RooArgSet* pdfDeps = pdf->getObservables(nset) ;
376  if (pdfDeps) {
377  supNSet.remove(*pdfDeps,kTRUE,kTRUE) ;
378  delete pdfDeps ;
379  }
380 
381  // Remove coef dependents
382  RooArgSet* coefDeps = coef ? coef->getObservables(nset) : 0 ;
383  if (coefDeps) {
384  supNSet.remove(*coefDeps,kTRUE,kTRUE) ;
385  delete coefDeps ;
386  }
387 
388  RooAbsReal* snorm ;
389  TString name(GetName()) ;
390  name.Append("_") ;
391  name.Append(pdf->GetName()) ;
392  name.Append("_SupNorm") ;
393  if (supNSet.getSize()>0) {
394  snorm = new RooRealIntegral(name,"Supplemental normalization integral",RooRealConstant::value(1.0),supNSet) ;
395  } else {
396  snorm = new RooRealVar(name,"Unit Supplemental normalization integral",1.0) ;
397  }
398  cache->_suppNormList.addOwned(*snorm) ;
399  }
400 
401  delete fullDepList ;
402 
403  if (_verboseEval>1) {
404  cxcoutD(Caching) << "RooAddModel::syncSuppNormList(" << GetName() << ") synching supplemental normalization list for norm" << (nset?*nset:RooArgSet()) << endl ;
405  if (dologD(Caching)) {
406  cache->_suppNormList.Print("v") ;
407  }
408  }
409 
410 
411  // *** PART 2 : Create projection coefficients ***
412 
413  // If no projections required stop here
414  if (!_projectCoefs || _basis!=0 ) {
415  _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
416  return cache ;
417  }
418 
419 
420  // Reduce iset/nset to actual dependents of this PDF
421  RooArgSet* nset2 = nset ? getObservables(nset) : new RooArgSet() ;
422 
423  // Check if requested transformation is not identity
424  if (!nset2->equals(_refCoefNorm) || _refCoefRangeName !=0 || rangeName !=0 ) {
425 
426  coutI(Caching) << "RooAddModel::syncCoefProjList(" << GetName() << ") creating coefficient projection integrals" << endl ;
427  ccoutI(Caching) << " from current normalization: " ; nset2->Print("1") ;
428  ccoutI(Caching) << " with current range: " << (rangeName?rangeName:"<none>") << endl ;
429  ccoutI(Caching) << " to reference normalization: " ; _refCoefNorm.Print("1") ;
430  ccoutI(Caching) << " with reference range: " << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"<none>") << endl ;
431 
432  // Recalculate projection integrals of PDFs
433  _pdfIter->Reset() ;
434  RooAbsPdf* thePdf ;
435 
436  while((thePdf=(RooAbsPdf*)_pdfIter->Next())) {
437 
438  // Calculate projection integral
439  RooAbsReal* pdfProj ;
440  if (!nset2->equals(_refCoefNorm)) {
441  pdfProj = thePdf->createIntegral(*nset2,_refCoefNorm) ;
442  pdfProj->setOperMode(operMode()) ;
443  } else {
444  TString name(GetName()) ;
445  name.Append("_") ;
446  name.Append(thePdf->GetName()) ;
447  name.Append("_ProjectNorm") ;
448  pdfProj = new RooRealVar(name,"Unit Projection normalization integral",1.0) ;
449  }
450 
451  cache->_projList.addOwned(*pdfProj) ;
452 
453  // Calculation optional supplemental normalization term
454  RooArgSet supNormSet(_refCoefNorm) ;
455  RooArgSet* deps = thePdf->getParameters(RooArgSet()) ;
456  supNormSet.remove(*deps,kTRUE,kTRUE) ;
457  delete deps ;
458 
459  RooAbsReal* snorm ;
460  TString name(GetName()) ;
461  name.Append("_") ;
462  name.Append(thePdf->GetName()) ;
463  name.Append("_ProjSupNorm") ;
464  if (supNormSet.getSize()>0) {
465  snorm = new RooRealIntegral(name,"Projection Supplemental normalization integral",
466  RooRealConstant::value(1.0),supNormSet) ;
467  } else {
468  snorm = new RooRealVar(name,"Unit Projection Supplemental normalization integral",1.0) ;
469  }
470  cache->_suppProjList.addOwned(*snorm) ;
471 
472  // Calculate reference range adjusted projection integral
473  RooAbsReal* rangeProj1 ;
474  if (_refCoefRangeName && _refCoefNorm.getSize()>0) {
475  rangeProj1 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,RooNameReg::str(_refCoefRangeName)) ;
476  rangeProj1->setOperMode(operMode()) ;
477  } else {
478  TString theName(GetName()) ;
479  theName.Append("_") ;
480  theName.Append(thePdf->GetName()) ;
481  theName.Append("_RangeNorm1") ;
482  rangeProj1 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
483  }
484  cache->_refRangeProjList.addOwned(*rangeProj1) ;
485 
486 
487  // Calculate range adjusted projection integral
488  RooAbsReal* rangeProj2 ;
489  if (rangeName && _refCoefNorm.getSize()>0) {
490  rangeProj2 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,rangeName) ;
491  rangeProj2->setOperMode(operMode()) ;
492  } else {
493  TString theName(GetName()) ;
494  theName.Append("_") ;
495  theName.Append(thePdf->GetName()) ;
496  theName.Append("_RangeNorm2") ;
497  rangeProj2 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
498  }
499  cache->_rangeProjList.addOwned(*rangeProj2) ;
500 
501  }
502 
503  }
504 
505  delete nset2 ;
506 
507  _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
508 
509  return cache ;
510 }
511 
512 
513 
514 ////////////////////////////////////////////////////////////////////////////////
515 /// Update the coefficient values in the given cache element: calculate new remainder
516 /// fraction, normalize fractions obtained from extended ML terms to unity and
517 /// multiply these the various range and dimensional corrections needed in the
518 /// current use context
519 
520 void RooAddModel::updateCoefficients(CacheElem& cache, const RooArgSet* nset) const
521 {
522  // cxcoutD(ChangeTracking) << "RooAddModel::updateCoefficients(" << GetName() << ") update coefficients" << endl ;
523 
524  Int_t i ;
525 
526  // Straight coefficients
527  if (_allExtendable) {
528 
529  // coef[i] = expectedEvents[i] / SUM(expectedEvents)
530  Double_t coefSum(0) ;
531  for (i=0 ; i<_pdfList.getSize() ; i++) {
532  _coefCache[i] = ((RooAbsPdf*)_pdfList.at(i))->expectedEvents(_refCoefNorm.getSize()>0?&_refCoefNorm:nset) ;
533  coefSum += _coefCache[i] ;
534  }
535  if (coefSum==0.) {
536  coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName() << ") WARNING: total number of expected events is 0" << endl ;
537  } else {
538  for (i=0 ; i<_pdfList.getSize() ; i++) {
539  _coefCache[i] /= coefSum ;
540  }
541  }
542 
543  } else {
544  if (_haveLastCoef) {
545 
546  // coef[i] = coef[i] / SUM(coef)
547  Double_t coefSum(0) ;
548  for (i=0 ; i<_coefList.getSize() ; i++) {
549  _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
550  coefSum += _coefCache[i] ;
551  }
552  for (i=0 ; i<_coefList.getSize() ; i++) {
553  _coefCache[i] /= coefSum ;
554  }
555  } else {
556 
557  // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
558  Double_t lastCoef(1) ;
559  for (i=0 ; i<_coefList.getSize() ; i++) {
560  _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
561  cxcoutD(Caching) << "SYNC: orig coef[" << i << "] = " << _coefCache[i] << endl ;
562  lastCoef -= _coefCache[i] ;
563  }
564  _coefCache[_coefList.getSize()] = lastCoef ;
565  cxcoutD(Caching) << "SYNC: orig coef[" << _coefList.getSize() << "] = " << _coefCache[_coefList.getSize()] << endl ;
566 
567 
568  // Warn about coefficient degeneration
569  if ((lastCoef<-1e-05 || (lastCoef-1)>1e-5) && _coefErrCount-->0) {
570  coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName()
571  << " WARNING: sum of PDF coefficients not in range [0-1], value="
572  << 1-lastCoef << endl ;
573  if (_coefErrCount==0) {
574  coutW(Eval) << " (no more will be printed)" << endl ;
575  }
576  }
577  }
578  }
579 
580 
581 
582  // Stop here if not projection is required or needed
583  if ((!_projectCoefs) || cache._projList.getSize()==0) {
584  // cout << "SYNC no projection required rangeName = " << (rangeName?rangeName:"<none>") << endl ;
585  return ;
586  }
587 
588  // Adjust coefficients for given projection
589  Double_t coefSum(0) ;
590  for (i=0 ; i<_pdfList.getSize() ; i++) {
591  GlobalSelectComponentRAII compRAII(true);
592 
593  RooAbsReal* pp = ((RooAbsReal*)cache._projList.at(i)) ;
594  RooAbsReal* sn = ((RooAbsReal*)cache._suppProjList.at(i)) ;
595  RooAbsReal* r1 = ((RooAbsReal*)cache._refRangeProjList.at(i)) ;
596  RooAbsReal* r2 = ((RooAbsReal*)cache._rangeProjList.at(i)) ;
597 
598  if (dologD(Eval)) {
599  cxcoutD(Eval) << "pp = " << pp->GetName() << endl
600  << "sn = " << sn->GetName() << endl
601  << "r1 = " << r1->GetName() << endl
602  << "r2 = " << r2->GetName() << endl ;
603  r1->printStream(ccoutD(Eval),kName|kArgs|kValue,kSingleLine) ;
604  r1->printCompactTree(ccoutD(Eval)) ;
605  }
606 
607  Double_t proj = pp->getVal()/sn->getVal()*(r2->getVal()/r1->getVal()) ;
608 
609  _coefCache[i] *= proj ;
610  coefSum += _coefCache[i] ;
611  }
612  for (i=0 ; i<_pdfList.getSize() ; i++) {
613  _coefCache[i] /= coefSum ;
614 // cout << "POST-SYNC coef[" << i << "] = " << _coefCache[i] << endl ;
615  }
616 
617 
618 
619 }
620 
621 
622 
623 ////////////////////////////////////////////////////////////////////////////////
624 /// Calculate the current value
625 
626 Double_t RooAddModel::evaluate() const
627 {
628  const RooArgSet* nset = _normSet ;
629  CacheElem* cache = getProjCache(nset) ;
630 
631  updateCoefficients(*cache,nset) ;
632 
633 
634  // Do running sum of coef/pdf pairs, calculate lastCoef.
635  _pdfIter->Reset() ;
636  _coefIter->Reset() ;
637  RooAbsPdf* pdf ;
638 
639  Double_t snormVal ;
640  Double_t value(0) ;
641  Int_t i(0) ;
642  while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
643  if (_coefCache[i]!=0.) {
644  snormVal = nset ? ((RooAbsReal*)cache->_suppNormList.at(i))->getVal() : 1.0 ;
645  Double_t pdfVal = pdf->getVal(nset) ;
646  // Double_t pdfNorm = pdf->getNorm(nset) ;
647  if (pdf->isSelectedComp()) {
648  value += pdfVal*_coefCache[i]/snormVal ;
649  cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
650  << pdf->GetName() << "] " << pdfVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
651  }
652  }
653  i++ ;
654  }
655 
656  return value ;
657 }
658 
659 
660 
661 ////////////////////////////////////////////////////////////////////////////////
662 /// Reset error counter to given value, limiting the number
663 /// of future error messages for this pdf to 'resetValue'
664 
665 void RooAddModel::resetErrorCounters(Int_t resetValue)
666 {
667  RooAbsPdf::resetErrorCounters(resetValue) ;
668  _coefErrCount = resetValue ;
669 }
670 
671 
672 
673 ////////////////////////////////////////////////////////////////////////////////
674 /// Check if PDF is valid for given normalization set.
675 /// Coeffient and PDF must be non-overlapping, but pdf-coefficient
676 /// pairs may overlap each other
677 
678 Bool_t RooAddModel::checkObservables(const RooArgSet* nset) const
679 {
680  Bool_t ret(kFALSE) ;
681 
682  _pdfIter->Reset() ;
683  _coefIter->Reset() ;
684  RooAbsReal* coef ;
685  RooAbsReal* pdf ;
686  while((coef=(RooAbsReal*)_coefIter->Next())) {
687  pdf = (RooAbsReal*)_pdfIter->Next() ;
688  if (pdf->observableOverlaps(nset,*coef)) {
689  coutE(InputArguments) << "RooAddModel::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName()
690  << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
691  ret = kTRUE ;
692  }
693  }
694 
695  return ret ;
696 }
697 
698 
699 
700 ////////////////////////////////////////////////////////////////////////////////
701 
702 Int_t RooAddModel::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars,
703  const RooArgSet* normSet, const char* rangeName) const
704 {
705  if (_forceNumInt) return 0 ;
706 
707  // Declare that we can analytically integrate all requested observables
708  analVars.add(allVars) ;
709 
710  // Retrieve (or create) the required component integral list
711  Int_t code ;
712  RooArgList *cilist ;
713  getCompIntList(normSet,&allVars,cilist,code,rangeName) ;
714 
715  return code+1 ;
716 
717 }
718 
719 
720 
721 ////////////////////////////////////////////////////////////////////////////////
722 /// Check if this configuration was created before
723 
724 void RooAddModel::getCompIntList(const RooArgSet* nset, const RooArgSet* iset, pRooArgList& compIntList, Int_t& code, const char* isetRangeName) const
725 {
726  Int_t sterileIdx(-1) ;
727 
728  IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObj(nset,iset,&sterileIdx,RooNameReg::ptr(isetRangeName)) ;
729  if (cache) {
730  code = _intCacheMgr.lastIndex() ;
731  compIntList = &cache->_intList ;
732 
733  return ;
734  }
735 
736  // Create containers for partial integral components to be generated
737  cache = new IntCacheElem ;
738 
739  // Fill Cache
740  _pdfIter->Reset() ;
741  RooResolutionModel* model ;
742  while ((model=(RooResolutionModel*)_pdfIter->Next())) {
743  RooAbsReal* intPdf = model->createIntegral(*iset,nset,0,isetRangeName) ;
744  cache->_intList.addOwned(*intPdf) ;
745  }
746 
747  // Store the partial integral list and return the assigned code ;
748  code = _intCacheMgr.setObj(nset,iset,(RooAbsCacheElement*)cache,RooNameReg::ptr(isetRangeName)) ;
749 
750  // Fill references to be returned
751  compIntList = &cache->_intList ;
752 }
753 
754 
755 
756 ////////////////////////////////////////////////////////////////////////////////
757 /// Return analytical integral defined by given scenario code
758 
759 Double_t RooAddModel::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
760 {
761  // No integration scenario
762  if (code==0) {
763  return getVal(normSet) ;
764  }
765 
766  // Partial integration scenarios
767  IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObjByIndex(code-1) ;
768 
769  RooArgList* compIntList ;
770 
771  // If cache has been sterilized, revive this slot
772  if (cache==0) {
773  RooArgSet* vars = getParameters(RooArgSet()) ;
774  RooArgSet* nset = _intCacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
775  RooArgSet* iset = _intCacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
776 
777  Int_t code2(-1) ;
778  getCompIntList(nset,iset,compIntList,code2,rangeName) ;
779 
780  delete vars ;
781  delete nset ;
782  delete iset ;
783  } else {
784 
785  compIntList = &cache->_intList ;
786 
787  }
788 
789  // Calculate the current value
790  const RooArgSet* nset = _normSet ;
791  CacheElem* pcache = getProjCache(nset) ;
792 
793  updateCoefficients(*pcache,nset) ;
794 
795  // Do running sum of coef/pdf pairs, calculate lastCoef.
796  TIterator* compIntIter = compIntList->createIterator() ;
797  _coefIter->Reset() ;
798  RooAbsReal* pdfInt ;
799 
800  Double_t snormVal ;
801  Double_t value(0) ;
802  Int_t i(0) ;
803  while((pdfInt = (RooAbsReal*)compIntIter->Next())) {
804  if (_coefCache[i]!=0.) {
805  snormVal = nset ? ((RooAbsReal*)pcache->_suppNormList.at(i))->getVal() : 1.0 ;
806  Double_t intVal = pdfInt->getVal(nset) ;
807  value += intVal*_coefCache[i]/snormVal ;
808  cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ") value += ["
809  << pdfInt->GetName() << "] " << intVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
810  }
811  i++ ;
812  }
813 
814  delete compIntIter ;
815 
816  return value ;
817 
818 }
819 
820 
821 
822 ////////////////////////////////////////////////////////////////////////////////
823 /// Return the number of expected events, which is either the sum of all coefficients
824 /// or the sum of the components extended terms
825 
826 Double_t RooAddModel::expectedEvents(const RooArgSet* nset) const
827 {
828  Double_t expectedTotal(0.0);
829  RooAbsPdf* pdf ;
830 
831  if (_allExtendable) {
832 
833  // Sum of the extended terms
834  _pdfIter->Reset() ;
835  while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
836  expectedTotal += pdf->expectedEvents(nset) ;
837  }
838 
839  } else {
840 
841  // Sum the coefficients
842  _coefIter->Reset() ;
843  RooAbsReal* coef ;
844  while((coef=(RooAbsReal*)_coefIter->Next())) {
845  expectedTotal += coef->getVal() ;
846  }
847  }
848 
849  return expectedTotal;
850 }
851 
852 
853 
854 ////////////////////////////////////////////////////////////////////////////////
855 /// Interface function used by test statistics to freeze choice of observables
856 /// for interpretation of fraction coefficients
857 
858 void RooAddModel::selectNormalization(const RooArgSet* depSet, Bool_t force)
859 {
860  if (!force && _refCoefNorm.getSize()!=0) {
861  return ;
862  }
863 
864  if (!depSet) {
865  fixCoefNormalization(RooArgSet()) ;
866  return ;
867  }
868 
869  RooArgSet* myDepSet = getObservables(depSet) ;
870  fixCoefNormalization(*myDepSet) ;
871  delete myDepSet ;
872 }
873 
874 
875 
876 ////////////////////////////////////////////////////////////////////////////////
877 /// Interface function used by test statistics to freeze choice of range
878 /// for interpretation of fraction coefficients
879 
880 void RooAddModel::selectNormalizationRange(const char* rangeName, Bool_t force)
881 {
882  if (!force && _refCoefRangeName) {
883  return ;
884  }
885 
886  fixCoefRange(rangeName) ;
887 }
888 
889 
890 
891 ////////////////////////////////////////////////////////////////////////////////
892 /// Return specialized context to efficiently generate toy events from RooAddPdfs
893 
894 RooAbsGenContext* RooAddModel::genContext(const RooArgSet &vars, const RooDataSet *prototype,
895  const RooArgSet* auxProto, Bool_t verbose) const
896 {
897  return new RooAddGenContext(*this,vars,prototype,auxProto,verbose) ;
898 }
899 
900 
901 
902 ////////////////////////////////////////////////////////////////////////////////
903 /// Direct generation is safe if all components say so
904 
905 Bool_t RooAddModel::isDirectGenSafe(const RooAbsArg& arg) const
906 {
907  _pdfIter->Reset() ;
908  RooAbsPdf* pdf ;
909  while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
910  if (!pdf->isDirectGenSafe(arg)) {
911  return kFALSE ;
912  }
913  }
914  return kTRUE ;
915 }
916 
917 
918 
919 ////////////////////////////////////////////////////////////////////////////////
920 /// Return pseud-code that indicates if all components can do internal generation (1) or not (0)
921 
922 Int_t RooAddModel::getGenerator(const RooArgSet& directVars, RooArgSet &/*generateVars*/, Bool_t /*staticInitOK*/) const
923 {
924  _pdfIter->Reset() ;
925  RooAbsPdf* pdf ;
926  while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
927  RooArgSet tmp ;
928  if (pdf->getGenerator(directVars,tmp)==0) {
929  return 0 ;
930  }
931  }
932  return 1 ;
933 }
934 
935 
936 
937 
938 ////////////////////////////////////////////////////////////////////////////////
939 /// This function should never be called as RooAddModel implements a custom generator context
940 
941 void RooAddModel::generateEvent(Int_t /*code*/)
942 {
943  assert(0) ;
944 }
945 
946 
947 
948 
949 ////////////////////////////////////////////////////////////////////////////////
950 /// List all RooAbsArg derived contents in this cache element
951 
952 RooArgList RooAddModel::CacheElem::containedArgs(Action)
953 {
954  RooArgList allNodes;
955  allNodes.add(_projList) ;
956  allNodes.add(_suppProjList) ;
957  allNodes.add(_refRangeProjList) ;
958  allNodes.add(_rangeProjList) ;
959 
960  return allNodes ;
961 }
962 
963 
964 
965 ////////////////////////////////////////////////////////////////////////////////
966 /// List all RooAbsArg derived contents in this cache element
967 
968 RooArgList RooAddModel::IntCacheElem::containedArgs(Action)
969 {
970  RooArgList allNodes(_intList) ;
971  return allNodes ;
972 }
973 
974 
975 ////////////////////////////////////////////////////////////////////////////////
976 /// Customized printing of arguments of a RooAddModel to more intuitively reflect the contents of the
977 /// product operator construction
978 
979 void RooAddModel::printMetaArgs(ostream& os) const
980 {
981  _pdfIter->Reset() ;
982  _coefIter->Reset() ;
983 
984  Bool_t first(kTRUE) ;
985 
986  os << "(" ;
987  RooAbsArg* coef, *pdf ;
988  while((coef=(RooAbsArg*)_coefIter->Next())) {
989  if (!first) {
990  os << " + " ;
991  } else {
992  first = kFALSE ;
993  }
994  pdf=(RooAbsArg*)_pdfIter->Next() ;
995  os << coef->GetName() << " * " << pdf->GetName() ;
996  }
997  pdf = (RooAbsArg*) _pdfIter->Next() ;
998  if (pdf) {
999  os << " + [%] * " << pdf->GetName() ;
1000  }
1001  os << ") " ;
1002 }
1003