RooAddPdf.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooAddPdf.cxx 36230 2010-10-09 20:21:02Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 // 
00019 // RooAddPdf is an efficient implementation of a sum of PDFs of the form 
00020 //
00021 //  c_1*PDF_1 + c_2*PDF_2 + ... c_n*PDF_n 
00022 //
00023 // or 
00024 //
00025 //  c_1*PDF_1 + c_2*PDF_2 + ... (1-sum(c_1...c_n-1))*PDF_n 
00026 //
00027 // The first form is for extended likelihood fits, where the
00028 // expected number of events is Sum(i) c_i. The coefficients c_i
00029 // can either be explicitly provided, or, if all components support
00030 // extended likelihood fits, they can be calculated the contribution
00031 // of each PDF to the total number of expected events.
00032 //
00033 // In the second form, the sum of the coefficients is enforced to be one,
00034 // and the coefficient of the last PDF is calculated from that condition.
00035 //
00036 // It is also possible to parameterize the coefficients recursively
00037 // 
00038 // c1*PDF_1 + (1-c1)(c2*PDF_2 + (1-c2)*(c3*PDF_3 + ....))
00039 //
00040 // In this form the sum of the coefficients is always less than 1.0
00041 // for all possible values of the individual coefficients between 0 and 1.
00042 // 
00043 // RooAddPdf relies on each component PDF to be normalized and will perform 
00044 // no normalization other than calculating the proper last coefficient c_n, if requested.
00045 // An (enforced) condition for this assuption is that each PDF_i is independent
00046 // of each coefficient_i.
00047 //
00048 // 
00049 
00050 #include "RooFit.h"
00051 #include "RooMsgService.h"
00052 
00053 #include "TIterator.h"
00054 #include "TIterator.h"
00055 #include "TList.h"
00056 #include "RooAddPdf.h"
00057 #include "RooDataSet.h"
00058 #include "RooRealProxy.h"
00059 #include "RooPlot.h"
00060 #include "RooRealVar.h"
00061 #include "RooAddGenContext.h"
00062 #include "RooRealConstant.h"
00063 #include "RooNameReg.h"
00064 #include "RooMsgService.h"
00065 #include "RooRecursiveFraction.h"
00066 #include "RooGlobalFunc.h"
00067 #include "RooRealIntegral.h"
00068 
00069 #include "Riostream.h"
00070 #include <algorithm>
00071 
00072 
00073 ClassImp(RooAddPdf)
00074 ;
00075 
00076 
00077 //_____________________________________________________________________________
00078 RooAddPdf::RooAddPdf() :
00079   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00080   _refCoefRangeName(0),
00081   _codeReg(10),
00082   _snormList(0)
00083 {
00084   // Default constructor used for persistence
00085 
00086   _pdfIter   = _pdfList.createIterator() ;
00087   _coefIter  = _coefList.createIterator() ;
00088 
00089   _coefCache = new Double_t[10] ;
00090   _coefErrCount = _errorCount ;
00091 }
00092 
00093 
00094 
00095 //_____________________________________________________________________________
00096 RooAddPdf::RooAddPdf(const char *name, const char *title) :
00097   RooAbsPdf(name,title), 
00098   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00099   _refCoefRangeName(0),
00100   _projectCoefs(kFALSE),
00101   _projCacheMgr(this,10),
00102   _codeReg(10),
00103   _pdfList("!pdfs","List of PDFs",this),
00104   _coefList("!coefficients","List of coefficients",this),
00105   _snormList(0),
00106   _haveLastCoef(kFALSE),
00107   _allExtendable(kFALSE)
00108 {
00109   // Dummy constructor 
00110   _pdfIter   = _pdfList.createIterator() ;
00111   _coefIter  = _coefList.createIterator() ;
00112 
00113   _coefCache = new Double_t[10] ;
00114   _coefErrCount = _errorCount ;
00115 
00116 }
00117 
00118 
00119 
00120 //_____________________________________________________________________________
00121 RooAddPdf::RooAddPdf(const char *name, const char *title,
00122                      RooAbsPdf& pdf1, RooAbsPdf& pdf2, RooAbsReal& coef1) : 
00123   RooAbsPdf(name,title),
00124   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00125   _refCoefRangeName(0),
00126   _projectCoefs(kFALSE),
00127   _projCacheMgr(this,10),
00128   _codeReg(10),
00129   _pdfList("!pdfs","List of PDFs",this),
00130   _coefList("!coefficients","List of coefficients",this),
00131   _haveLastCoef(kFALSE),
00132   _allExtendable(kFALSE)
00133 {
00134   // Constructor with two PDFs and one coefficient
00135 
00136   _pdfIter  = _pdfList.createIterator() ;
00137   _coefIter = _coefList.createIterator() ;
00138 
00139   _pdfList.add(pdf1) ;  
00140   _pdfList.add(pdf2) ;
00141   _coefList.add(coef1) ;
00142 
00143   _coefCache = new Double_t[_pdfList.getSize()] ;
00144   _coefErrCount = _errorCount ;
00145 }
00146 
00147 
00148 
00149 //_____________________________________________________________________________
00150 RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, Bool_t recursiveFractions) :
00151   RooAbsPdf(name,title),
00152   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00153   _refCoefRangeName(0),
00154   _projectCoefs(kFALSE),
00155   _projCacheMgr(this,10),
00156   _codeReg(10),
00157   _pdfList("!pdfs","List of PDFs",this),
00158   _coefList("!coefficients","List of coefficients",this),
00159   _haveLastCoef(kFALSE),
00160   _allExtendable(kFALSE)
00161 { 
00162   // Generic constructor from list of PDFs and list of coefficients.
00163   // Each pdf list element (i) is paired with coefficient list element (i).
00164   // The number of coefficients must be either equal to the number of PDFs,
00165   // in which case extended MLL fitting is enabled, or be one less.
00166   //
00167   // All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal
00168   //
00169   // If the recursiveFraction flag is true, the coefficients are interpreted as recursive
00170   // coefficients as explained in the class description.
00171 
00172   if (inPdfList.getSize()>inCoefList.getSize()+1 || inPdfList.getSize()<inCoefList.getSize()) {
00173     coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() 
00174                           << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
00175     assert(0) ;
00176   }
00177 
00178   if (recursiveFractions && inPdfList.getSize()!=inCoefList.getSize()+1) {
00179     coutW(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() 
00180                           << ") WARNING inconsistent input: recursive fractions options can only be used if Npdf=Ncoef+1, ignoring recursive fraction setting" << endl ;
00181   }
00182 
00183 
00184   _pdfIter  = _pdfList.createIterator() ;
00185   _coefIter = _coefList.createIterator() ;
00186  
00187   // Constructor with N PDFs and N or N-1 coefs
00188   TIterator* pdfIter = inPdfList.createIterator() ;
00189   TIterator* coefIter = inCoefList.createIterator() ;
00190   RooAbsPdf* pdf ;
00191   RooAbsReal* coef ;
00192 
00193   RooArgList partinCoefList ;
00194 
00195   Bool_t first(kTRUE) ;
00196 
00197   while((coef = (RooAbsPdf*)coefIter->Next())) {
00198     pdf = (RooAbsPdf*) pdfIter->Next() ;
00199     if (!pdf) {
00200       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() 
00201                             << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
00202       assert(0) ;
00203     }
00204     if (!dynamic_cast<RooAbsReal*>(coef)) {
00205       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal, ignored" << endl ;
00206       continue ;
00207     }
00208     if (!dynamic_cast<RooAbsReal*>(pdf)) {
00209       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName() << " is not of type RooAbsPdf, ignored" << endl ;
00210       continue ;
00211     }
00212     _pdfList.add(*pdf) ;
00213 
00214     // Process recursive fraction mode separately
00215     if (recursiveFractions) {
00216       partinCoefList.add(*coef) ;
00217       if (first) {      
00218 
00219         // The first fraction is the first plain fraction
00220         first = kFALSE ;
00221         _coefList.add(*coef) ;
00222 
00223       } else {
00224 
00225         // The i-th recursive fraction = (1-f1)*(1-f2)*...(fi) and is calculated from the list (f1,...,fi) by RooRecursiveFraction)
00226         RooAbsReal* rfrac = new RooRecursiveFraction(Form("%s_recursive_fraction_%s",GetName(),pdf->GetName()),"Recursive Fraction",partinCoefList) ;
00227         addOwnedComponents(*rfrac) ;
00228         _coefList.add(*rfrac) ;         
00229 
00230       }
00231 
00232     } else {
00233       _coefList.add(*coef) ;    
00234     }
00235   }
00236 
00237   pdf = (RooAbsPdf*) pdfIter->Next() ;
00238   if (pdf) {
00239     if (!dynamic_cast<RooAbsReal*>(pdf)) {
00240       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") last pdf " << coef->GetName() << " is not of type RooAbsPdf, fatal error" << endl ;
00241       assert(0) ;
00242     }
00243     _pdfList.add(*pdf) ;  
00244 
00245     // Process recursive fractions mode
00246     if (recursiveFractions) {
00247 
00248       // The last recursive fraction = (1-f1)*(1-f2)*...(1-fN) and is calculated from the list (f1,...,fN,1) by RooRecursiveFraction)
00249       partinCoefList.add(RooFit::RooConst(1)) ;
00250       RooAbsReal* rfrac = new RooRecursiveFraction(Form("%s_recursive_fraction_%s",GetName(),pdf->GetName()),"Recursive Fraction",partinCoefList) ;
00251       addOwnedComponents(*rfrac) ;
00252       _coefList.add(*rfrac) ;           
00253 
00254       // In recursive mode we always have Ncoef=Npdf
00255       _haveLastCoef=kTRUE ;
00256     }
00257 
00258   } else {
00259     _haveLastCoef=kTRUE ;
00260   }
00261 
00262   delete pdfIter ;
00263   delete coefIter  ;
00264 
00265   _coefCache = new Double_t[_pdfList.getSize()] ;
00266   _coefErrCount = _errorCount ;
00267 }
00268 
00269 
00270 
00271 //_____________________________________________________________________________
00272 RooAddPdf::RooAddPdf(const char *name, const char *title, const RooArgList& inPdfList) :
00273   RooAbsPdf(name,title),
00274   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00275   _refCoefRangeName(0),
00276   _projectCoefs(kFALSE),
00277   _projCacheMgr(this,10),
00278   _pdfList("!pdfs","List of PDFs",this),
00279   _coefList("!coefficients","List of coefficients",this),
00280   _haveLastCoef(kFALSE),
00281   _allExtendable(kTRUE)
00282 { 
00283   // Generic constructor from list of extended PDFs. There are no coefficients as the expected
00284   // number of events from each components determine the relative weight of the PDFs.
00285   // 
00286   // All PDFs must inherit from RooAbsPdf. 
00287 
00288   _pdfIter  = _pdfList.createIterator() ;
00289   _coefIter = _coefList.createIterator() ;
00290  
00291   // Constructor with N PDFs 
00292   TIterator* pdfIter = inPdfList.createIterator() ;
00293   RooAbsPdf* pdf ;
00294   while((pdf = (RooAbsPdf*) pdfIter->Next())) {
00295     
00296     if (!dynamic_cast<RooAbsReal*>(pdf)) {
00297       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName() << " is not of type RooAbsPdf, ignored" << endl ;
00298       continue ;
00299     }
00300     if (!pdf->canBeExtended()) {
00301       coutE(InputArguments) << "RooAddPdf::RooAddPdf(" << GetName() << ") pdf " << pdf->GetName() << " is not extendable, ignored" << endl ;
00302       continue ;
00303     }
00304     _pdfList.add(*pdf) ;    
00305   }
00306 
00307   delete pdfIter ;
00308 
00309   _coefCache = new Double_t[_pdfList.getSize()] ;
00310   _coefErrCount = _errorCount ;
00311 }
00312 
00313 
00314 
00315 
00316 //_____________________________________________________________________________
00317 RooAddPdf::RooAddPdf(const RooAddPdf& other, const char* name) :
00318   RooAbsPdf(other,name),
00319   _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
00320   _refCoefRangeName((TNamed*)other._refCoefRangeName),
00321   _projectCoefs(other._projectCoefs),
00322   _projCacheMgr(other._projCacheMgr,this),
00323   _codeReg(other._codeReg),
00324   _pdfList("!pdfs",this,other._pdfList),
00325   _coefList("!coefficients",this,other._coefList),
00326   _haveLastCoef(other._haveLastCoef),
00327   _allExtendable(other._allExtendable)
00328 {
00329   // Copy constructor
00330 
00331   _pdfIter  = _pdfList.createIterator() ;
00332   _coefIter = _coefList.createIterator() ;
00333   _coefCache = new Double_t[_pdfList.getSize()] ;
00334   _coefErrCount = _errorCount ;
00335 }
00336 
00337 
00338 
00339 //_____________________________________________________________________________
00340 RooAddPdf::~RooAddPdf()
00341 {
00342   // Destructor
00343 
00344   delete _pdfIter ;
00345   delete _coefIter ;
00346 
00347   if (_coefCache) delete[] _coefCache ;
00348 }
00349 
00350 
00351 
00352 //_____________________________________________________________________________
00353 void RooAddPdf::fixCoefNormalization(const RooArgSet& refCoefNorm) 
00354 {
00355   // By default the interpretation of the fraction coefficients is
00356   // performed in the contextual choice of observables. This makes the
00357   // shape of the p.d.f explicitly dependent on the choice of
00358   // observables. This method instructs RooAddPdf to freeze the
00359   // interpretation of the coefficients to be done in the given set of
00360   // observables. If frozen, fractions are automatically transformed
00361   // from the reference normalization set to the contextual normalization
00362   // set by ratios of integrals
00363 
00364   if (refCoefNorm.getSize()==0) {
00365     _projectCoefs = kFALSE ;
00366     return ;
00367   }
00368   _projectCoefs = kTRUE ;  
00369 
00370   _refCoefNorm.removeAll() ;
00371   _refCoefNorm.add(refCoefNorm) ;
00372 
00373   _projCacheMgr.reset() ;
00374 }
00375 
00376 
00377 
00378 //_____________________________________________________________________________
00379 void RooAddPdf::fixCoefRange(const char* rangeName)
00380 {
00381   // By default the interpretation of the fraction coefficients is
00382   // performed in the default range. This make the shape of a RooAddPdf
00383   // explicitly dependent on the range of the observables. To allow
00384   // a range independent definition of the fraction this function
00385   // instructs RooAddPdf to freeze its interpretation in the given
00386   // named range. If the current normalization range is different
00387   // from the reference range, the appropriate fraction coefficients
00388   // are automically calculation from the reference fractions using
00389   // ratios if integrals
00390 
00391   _refCoefRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
00392   if (_refCoefRangeName) _projectCoefs = kTRUE ;
00393 }
00394 
00395 
00396 
00397 //_____________________________________________________________________________
00398 RooAddPdf::CacheElem* RooAddPdf::getProjCache(const RooArgSet* nset, const RooArgSet* iset, const char* rangeName) const
00399 {
00400   // Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated over iset
00401   // in range 'rangeName'. If cache element does not exist, create and fill it on the fly. The cache contains
00402   // suplemental normalization terms (in case not all added p.d.f.s have the same observables), projection
00403   // integrals to calculated transformed fraction coefficients when a frozen reference frame is provided
00404   // and projection integrals for similar transformations when a frozen reference range is provided.
00405 
00406 
00407   // Check if cache already exists 
00408   CacheElem* cache = (CacheElem*) _projCacheMgr.getObj(nset,iset,0,RooNameReg::ptr(rangeName)) ;
00409   if (cache) {
00410     return cache ;
00411   }
00412 
00413   //Create new cache 
00414   cache = new CacheElem ;
00415 
00416   // *** PART 1 : Create supplemental normalization list ***
00417 
00418   // Retrieve the combined set of dependents of this PDF ;
00419   RooArgSet *fullDepList = getObservables(nset) ;
00420   if (iset) {
00421     fullDepList->remove(*iset,kTRUE,kTRUE) ;
00422   }    
00423 
00424   // Fill with dummy unit RRVs for now
00425   _pdfIter->Reset() ;
00426   _coefIter->Reset() ;
00427   RooAbsPdf* pdf ;
00428   RooAbsReal* coef ;
00429   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {    
00430     coef=(RooAbsPdf*)_coefIter->Next() ;
00431 
00432     // Start with full list of dependents
00433     RooArgSet supNSet(*fullDepList) ;
00434 
00435     // Remove PDF dependents
00436     RooArgSet* pdfDeps = pdf->getObservables(nset) ;
00437     if (pdfDeps) {
00438       supNSet.remove(*pdfDeps,kTRUE,kTRUE) ;
00439       delete pdfDeps ; 
00440     }
00441 
00442     // Remove coef dependents
00443     RooArgSet* coefDeps = coef ? coef->getObservables(nset) : 0 ;
00444     if (coefDeps) {
00445       supNSet.remove(*coefDeps,kTRUE,kTRUE) ;
00446       delete coefDeps ;
00447     }
00448     
00449     RooAbsReal* snorm ;
00450     TString name(GetName()) ;
00451     name.Append("_") ;
00452     name.Append(pdf->GetName()) ;
00453     name.Append("_SupNorm") ;
00454     if (supNSet.getSize()>0) {
00455       snorm = new RooRealIntegral(name,"Supplemental normalization integral",RooRealConstant::value(1.0),supNSet) ;
00456       cxcoutD(Caching) << "RooAddPdf " << GetName() << " making supplemental normalization set " << supNSet << " for pdf component " << pdf->GetName() << endl ;
00457     } else {
00458       snorm = new RooRealVar(name,"Unit Supplemental normalization integral",1.0) ;
00459     }
00460     cache->_suppNormList.addOwned(*snorm) ;
00461   }
00462 
00463   delete fullDepList ;
00464     
00465   if (_verboseEval>1) {
00466     cxcoutD(Caching) << "RooAddPdf::syncSuppNormList(" << GetName() << ") synching supplemental normalization list for norm" << (nset?*nset:RooArgSet()) << endl ;
00467     if dologD(Caching) {
00468       cache->_suppNormList.Print("v") ;
00469     }
00470   }
00471 
00472 
00473   // *** PART 2 : Create projection coefficients ***
00474 
00475 //   cout << " this = " << this << " (" << GetName() << ")" << endl ;
00476 //   cout << "projectCoefs = " << (_projectCoefs?"T":"F") << endl ;
00477 //   cout << "_normRange.Length() = " << _normRange.Length() << endl ;
00478 
00479   // If no projections required stop here
00480   if (!_projectCoefs && !rangeName) {
00481     _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
00482 //     cout << " no projection required" << endl ;
00483     return cache ;
00484   }
00485 
00486 
00487 //   cout << "calculating projection" << endl ;
00488 
00489   // Reduce iset/nset to actual dependents of this PDF
00490   RooArgSet* nset2 = nset ? getObservables(nset) : new RooArgSet() ;
00491   cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset = " << (nset?*nset:RooArgSet()) << " nset2 = " << *nset2 << endl ;
00492 
00493   if (nset2->getSize()==0 && _refCoefNorm.getSize()!=0) {
00494     //cout << "WVE: evaluating RooAddPdf without normalization, but have reference normalization for coefficient definition" << endl ;
00495     
00496     nset2->add(_refCoefNorm) ;
00497     if (_refCoefRangeName) {
00498       rangeName = RooNameReg::str(_refCoefRangeName) ;
00499     }
00500   }
00501 
00502 
00503   // Check if requested transformation is not identity 
00504   if (!nset2->equals(_refCoefNorm) || _refCoefRangeName !=0 || rangeName !=0 || _normRange.Length()>0) {
00505    
00506     cxcoutD(Caching) << "ALEX:     RooAddPdf::syncCoefProjList(" << GetName() << ") projecting coefficients from "
00507                    << *nset2 << (rangeName?":":"") << (rangeName?rangeName:"") 
00508                    << " to "  << ((_refCoefNorm.getSize()>0)?_refCoefNorm:*nset2) << (_refCoefRangeName?":":"") << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"") << endl ;
00509     
00510     // Recalculate projection integrals of PDFs 
00511     _pdfIter->Reset() ;
00512     RooAbsPdf* thePdf ;
00513 
00514     while((thePdf=(RooAbsPdf*)_pdfIter->Next())) {
00515 
00516       // Calculate projection integral
00517       RooAbsReal* pdfProj ;
00518       if (!nset2->equals(_refCoefNorm)) {
00519         pdfProj = thePdf->createIntegral(*nset2,_refCoefNorm,_normRange.Length()>0?_normRange.Data():0) ;
00520         pdfProj->setOperMode(operMode()) ;
00521         cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset2(" << *nset2 << ")!=_refCoefNorm(" << _refCoefNorm << ") --> pdfProj = " << pdfProj->GetName() << endl ;
00522       } else {
00523         TString name(GetName()) ;
00524         name.Append("_") ;
00525         name.Append(thePdf->GetName()) ;
00526         name.Append("_ProjectNorm") ;
00527         pdfProj = new RooRealVar(name,"Unit Projection normalization integral",1.0) ;
00528         cxcoutD(Caching) << "RooAddPdf(" << GetName() << ")::getPC nset2(" << *nset2 << ")==_refCoefNorm(" << _refCoefNorm << ") --> pdfProj = " << pdfProj->GetName() << endl ;
00529       }
00530 
00531       cache->_projList.addOwned(*pdfProj) ;
00532       cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") PP = " << pdfProj->GetName() << endl ;
00533 
00534       // Calculation optional supplemental normalization term
00535       RooArgSet supNormSet(_refCoefNorm) ;
00536       RooArgSet* deps = thePdf->getParameters(RooArgSet()) ;
00537       supNormSet.remove(*deps,kTRUE,kTRUE) ;
00538       delete deps ;
00539 
00540       RooAbsReal* snorm ;
00541       TString name(GetName()) ;
00542       name.Append("_") ;
00543       name.Append(thePdf->GetName()) ;
00544       name.Append("_ProjSupNorm") ;
00545       if (supNormSet.getSize()>0 && !nset2->equals(_refCoefNorm) ) {
00546         snorm = new RooRealIntegral(name,"Projection Supplemental normalization integral",
00547                                     RooRealConstant::value(1.0),supNormSet) ;
00548       } else {
00549         snorm = new RooRealVar(name,"Unit Projection Supplemental normalization integral",1.0) ;
00550       }
00551       cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") SN = " << snorm->GetName() << endl ;
00552       cache->_suppProjList.addOwned(*snorm) ;
00553 
00554       // Calculate reference range adjusted projection integral
00555       RooAbsReal* rangeProj1 ;
00556 
00557    //    cout << "ALEX >>>> RooAddPdf(" << GetName() << ")::getPC _refCoefRangeName WVE = "        
00558 //         <<(_refCoefRangeName?":":"") << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"")
00559 //         <<" _refCoefRangeName AK = "  << (_refCoefRangeName?_refCoefRangeName->GetName():"")            
00560 //         << " && _refCoefNorm" << _refCoefNorm << " with size = _refCoefNorm.getSize() " << _refCoefNorm.getSize() << endl ;
00561 
00562       if (_refCoefRangeName && _refCoefNorm.getSize()>0) {
00563         RooArgSet* tmp = thePdf->getObservables(_refCoefNorm) ;
00564         rangeProj1 = thePdf->createIntegral(*tmp,*tmp,RooNameReg::str(_refCoefRangeName)) ;
00565         //rangeProj1->setOperMode(operMode()) ;
00566         delete tmp ;
00567       } else {
00568         TString theName(GetName()) ;
00569         theName.Append("_") ;
00570         theName.Append(thePdf->GetName()) ;
00571         theName.Append("_RangeNorm1") ;
00572         rangeProj1 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
00573       }
00574       cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") R1 = " << rangeProj1->GetName() << endl ;
00575       cache->_refRangeProjList.addOwned(*rangeProj1) ;
00576       
00577 
00578       // Calculate range adjusted projection integral
00579       RooAbsReal* rangeProj2 ;
00580       cxcoutD(Caching) << "RooAddPdf::syncCoefProjList(" << GetName() << ") rangename = " << (rangeName?rangeName:"<null>") 
00581                        << " nset = " << (nset?*nset:RooArgSet()) << endl ;
00582       if (rangeName && _refCoefNorm.getSize()>0) {
00583 
00584         rangeProj2 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,rangeName) ;
00585         //rangeProj2->setOperMode(operMode()) ;
00586 
00587       } else if (_normRange.Length()>0) {
00588 
00589         RooArgSet* tmp = thePdf->getObservables(_refCoefNorm) ;
00590         rangeProj2 = thePdf->createIntegral(*tmp,*tmp,_normRange.Data()) ;
00591         delete tmp ;
00592 
00593       } else {
00594 
00595         TString theName(GetName()) ;
00596         theName.Append("_") ;
00597         theName.Append(thePdf->GetName()) ;
00598         theName.Append("_RangeNorm2") ;
00599         rangeProj2 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
00600 
00601       }
00602       cxcoutD(Caching) << " RooAddPdf::syncCoefProjList(" << GetName() << ") R2 = " << rangeProj2->GetName() << endl ;
00603       cache->_rangeProjList.addOwned(*rangeProj2) ;
00604 
00605     }               
00606 
00607   }
00608 
00609   delete nset2 ;
00610 
00611   _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
00612 
00613   return cache ;
00614 }
00615 
00616 
00617 //_____________________________________________________________________________
00618 void RooAddPdf::updateCoefficients(CacheElem& cache, const RooArgSet* nset) const 
00619 {
00620   // Update the coefficient values in the given cache element: calculate new remainder
00621   // fraction, normalize fractions obtained from extended ML terms to unity and
00622   // multiply these the various range and dimensional corrections needed in the
00623   // current use context
00624 
00625   // cxcoutD(ChangeTracking) << "RooAddPdf::updateCoefficients(" << GetName() << ") update coefficients" << endl ;
00626   
00627   Int_t i ;
00628 
00629   // Straight coefficients
00630   if (_allExtendable) {
00631     
00632     // coef[i] = expectedEvents[i] / SUM(expectedEvents)
00633     Double_t coefSum(0) ;
00634     for (i=0 ; i<_pdfList.getSize() ; i++) {
00635       _coefCache[i] = ((RooAbsPdf*)_pdfList.at(i))->expectedEvents(_refCoefNorm.getSize()>0?&_refCoefNorm:nset) ;
00636       coefSum += _coefCache[i] ;
00637     }
00638     if (coefSum==0.) {
00639       coutW(Eval) << "RooAddPdf::updateCoefCache(" << GetName() << ") WARNING: total number of expected events is 0" << endl ;
00640     } else {
00641       for (i=0 ; i<_pdfList.getSize() ; i++) {
00642         _coefCache[i] /= coefSum ;
00643       }                             
00644     }
00645     
00646   } else {
00647     if (_haveLastCoef) {
00648       
00649       // coef[i] = coef[i] / SUM(coef)
00650       Double_t coefSum(0) ;
00651       for (i=0 ; i<_coefList.getSize() ; i++) {
00652         _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
00653         coefSum += _coefCache[i] ;
00654       }         
00655       if (coefSum==0.) {
00656         coutW(Eval) << "RooAddPdf::updateCoefCache(" << GetName() << ") WARNING: sum of coefficients is zero 0" << endl ;
00657       } else {  
00658         for (i=0 ; i<_coefList.getSize() ; i++) {
00659           _coefCache[i] /= coefSum ;
00660         }                       
00661       }
00662     } else {
00663       
00664       // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
00665       Double_t lastCoef(1) ;
00666       for (i=0 ; i<_coefList.getSize() ; i++) {
00667         _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
00668         cxcoutD(Caching) << "SYNC: orig coef[" << i << "] = " << _coefCache[i] << endl ;
00669         lastCoef -= _coefCache[i] ;
00670       }                 
00671       _coefCache[_coefList.getSize()] = lastCoef ;
00672       cxcoutD(Caching) << "SYNC: orig coef[" << _coefList.getSize() << "] = " << _coefCache[_coefList.getSize()] << endl ;
00673       
00674       
00675       // Warn about coefficient degeneration
00676       if ((lastCoef<-1e-05 || (lastCoef-1)>1e-5) && _coefErrCount-->0) {
00677         coutW(Eval) << "RooAddPdf::updateCoefCache(" << GetName() 
00678                     << " WARNING: sum of PDF coefficients not in range [0-1], value=" 
00679                     << 1-lastCoef ; 
00680         if (_coefErrCount==0) {
00681           coutW(Eval) << " (no more will be printed)"  ;
00682         }
00683         coutW(Eval) << endl ;
00684       } 
00685     }
00686   }
00687 
00688   
00689 
00690 //    cout << "XXXX" << GetName() << "updateCoefs _projectCoefs = " << (_projectCoefs?"T":"F") << " cache._projList.getSize()= " << cache._projList.getSize() << endl ;
00691 
00692   // Stop here if not projection is required or needed
00693   if ((!_projectCoefs && _normRange.Length()==0) || cache._projList.getSize()==0) {
00694     //if (cache._projList.getSize()==0) {
00695 //     cout << GetName() << " SYNC no projection required rangeName = " << (_normRange.Length()>0?_normRange.Data():"<none>") << endl ;
00696     return ;
00697   }
00698 
00699 //    cout << "XXXX" << GetName() << " updateCoefs, applying correction" << endl ;
00700 //    cout << "PROJLIST = " << endl ;
00701 //    cache._projList.Print("v") ;
00702    
00703    
00704   // Adjust coefficients for given projection
00705   Double_t coefSum(0) ;
00706   for (i=0 ; i<_pdfList.getSize() ; i++) {
00707     Bool_t _tmp = _globalSelectComp ;
00708     RooAbsPdf::globalSelectComp(kTRUE) ;    
00709 
00710     RooAbsReal* pp = ((RooAbsReal*)cache._projList.at(i)) ; 
00711     RooAbsReal* sn = ((RooAbsReal*)cache._suppProjList.at(i)) ; 
00712     RooAbsReal* r1 = ((RooAbsReal*)cache._refRangeProjList.at(i)) ;
00713     RooAbsReal* r2 = ((RooAbsReal*)cache._rangeProjList.at(i)) ;
00714 
00715     Double_t proj = pp->getVal()/sn->getVal()*(r2->getVal()/r1->getVal()) ;  
00716     
00717     cxcoutD(Caching) << "ALEX:    RooAddPdf::updateCoef(" << GetName() << ") with nset = " << (nset?*nset:RooArgSet()) << "for pdf component #" << i << " = " << _pdfList.at(i)->GetName() << endl
00718          << "ALEX:   pp = " << pp->GetName() << " = " << pp->getVal() << endl 
00719          << "ALEX:   sn = " << sn->GetName() << " = " << sn->getVal() <<  endl 
00720          << "ALEX:   r1 = " << r1->GetName() << " = " << r1->getVal() <<  endl 
00721          << "ALEX:   r2 = " << r2->GetName() << " = " << r2->getVal() <<  endl 
00722          << "ALEX: proj = (" << pp->getVal() << "/" << sn->getVal() << ")*(" << r2->getVal() << "/" << r1->getVal() << ") = " << proj << endl ;
00723     
00724     
00725     
00726     RooAbsPdf::globalSelectComp(_tmp) ;
00727 
00728     _coefCache[i] *= proj ;
00729     coefSum += _coefCache[i] ;
00730   }
00731 
00732   for (i=0 ; i<_pdfList.getSize() ; i++) {
00733     _coefCache[i] /= coefSum ;
00734 //    _coefCache[i] *= rfrac ;
00735 //       cout << "POST-SYNC coef[" << i << "] = " << _coefCache[i] << endl ;
00736      cxcoutD(Caching) << " ALEX:   POST-SYNC coef[" << i << "] = " << _coefCache[i]  
00737          << " ( _coefCache[i]/coefSum = " << _coefCache[i]*coefSum << "/" << coefSum << " ) "<< endl ;
00738   }
00739    
00740 
00741   
00742 }
00743 
00744 
00745 
00746 //_____________________________________________________________________________
00747 Double_t RooAddPdf::evaluate() const 
00748 {
00749   // Calculate and return the current value
00750 
00751   const RooArgSet* nset = _normSet ; 
00752   cxcoutD(Caching) << "RooAddPdf::evaluate(" << GetName() << ") calling getProjCache with nset = " << nset << " = " << (nset?*nset:RooArgSet()) << endl ;
00753 
00754   if (nset==0 || nset->getSize()==0) {
00755     if (_refCoefNorm.getSize()!=0) {
00756       //cout << "RooAddPdf::eval(" << GetName() << ") WARNING NSET IS NULL, but have reference normalization, using that" << endl ;
00757       nset = &_refCoefNorm ;
00758     }
00759   }
00760 
00761   CacheElem* cache = getProjCache(nset) ;
00762   updateCoefficients(*cache,nset) ;
00763   
00764   // Do running sum of coef/pdf pairs, calculate lastCoef.
00765   _pdfIter->Reset() ;
00766   _coefIter->Reset() ;
00767   RooAbsPdf* pdf ;
00768 
00769   Double_t snormVal ;
00770   Double_t value(0) ;
00771   Int_t i(0) ;
00772   while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
00773     if (_coefCache[i]!=0.) {
00774       snormVal = nset ? ((RooAbsReal*)cache->_suppNormList.at(i))->getVal() : 1.0 ;
00775       Double_t pdfVal = pdf->getVal(nset) ;
00776       // Double_t pdfNorm = pdf->getNorm(nset) ;
00777       if (pdf->isSelectedComp()) {
00778         value += pdfVal*_coefCache[i]/snormVal ;
00779 //      cout << "RooAddPdf::EVALUATE(" << GetName() << ") adding pdf " << pdf->GetName() << " value = " << pdfVal << " coef = " << _coefCache[i] << endl ;
00780       }
00781     }
00782     i++ ;
00783   }
00784 
00785   return value ;
00786 }
00787 
00788 
00789 //_____________________________________________________________________________
00790 void RooAddPdf::resetErrorCounters(Int_t resetValue)
00791 {
00792   // Reset error counter to given value, limiting the number
00793   // of future error messages for this pdf to 'resetValue'
00794 
00795   RooAbsPdf::resetErrorCounters(resetValue) ;
00796   _coefErrCount = resetValue ;
00797 }
00798 
00799 
00800 
00801 //_____________________________________________________________________________
00802 Bool_t RooAddPdf::checkObservables(const RooArgSet* nset) const 
00803 {
00804   // Check if PDF is valid for given normalization set.
00805   // Coeffient and PDF must be non-overlapping, but pdf-coefficient 
00806   // pairs may overlap each other
00807 
00808   Bool_t ret(kFALSE) ;
00809 
00810   _pdfIter->Reset() ;
00811   _coefIter->Reset() ;
00812   RooAbsReal* coef ;
00813   RooAbsReal* pdf ;
00814   while((coef=(RooAbsReal*)_coefIter->Next())) {
00815     pdf = (RooAbsReal*)_pdfIter->Next() ;
00816     if (pdf->observableOverlaps(nset,*coef)) {
00817       coutE(InputArguments) << "RooAddPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName() 
00818                             << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
00819       ret = kTRUE ;
00820     }
00821   }
00822   
00823   return ret ;
00824 }
00825 
00826 
00827 //_____________________________________________________________________________
00828 Int_t RooAddPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
00829                                          const RooArgSet* normSet, const char* rangeName) const 
00830 {
00831   // Determine which part (if any) of given integral can be performed analytically.
00832   // If any analytical integration is possible, return integration scenario code
00833   //
00834   // RooAddPdf queries each component PDF for its analytical integration capability of the requested
00835   // set ('allVars'). It finds the largest common set of variables that can be integrated
00836   // by all components. If such a set exists, it reconfirms that each component is capable of
00837   // analytically integrating the common set, and combines the components individual integration
00838   // codes into a single integration code valid for RooAddPdf.
00839 
00840 
00841   RooArgSet* allDepVars = getObservables(allVars) ;
00842   RooArgSet allAnalVars(*allDepVars) ;
00843   delete allDepVars ;
00844 
00845   TIterator* avIter = allVars.createIterator() ;
00846 
00847   Int_t n(0) ;
00848 
00849   // First iteration, determine what each component can integrate analytically
00850   _pdfIter->Reset() ;
00851   RooAbsPdf* pdf ;
00852   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
00853     RooArgSet subAnalVars ;
00854     pdf->getAnalyticalIntegralWN(allVars,subAnalVars,normSet,rangeName) ;
00855 
00856     // Observables that cannot be integrated analytically by this component are dropped from the common list
00857     avIter->Reset() ;
00858     RooAbsArg* arg ;
00859     while((arg=(RooAbsArg*)avIter->Next())) {
00860       if (!subAnalVars.find(arg->GetName()) && pdf->dependsOn(*arg)) {
00861         allAnalVars.remove(*arg,kTRUE,kTRUE) ;
00862       } 
00863     }
00864     n++ ;
00865   }
00866 
00867   // If no observables can be integrated analytically, return code 0 here
00868   if (allAnalVars.getSize()==0) {
00869     delete avIter ;
00870     return 0 ;
00871   }
00872 
00873 
00874   // Now retrieve codes for integration over common set of analytically integrable observables for each component
00875   _pdfIter->Reset() ;
00876   n=0 ;
00877   Int_t* subCode = new Int_t[_pdfList.getSize()] ;
00878   Bool_t allOK(kTRUE) ;
00879   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
00880     RooArgSet subAnalVars ;
00881     RooArgSet* allAnalVars2 = pdf->getObservables(allAnalVars) ;
00882     subCode[n] = pdf->getAnalyticalIntegralWN(*allAnalVars2,subAnalVars,normSet,rangeName) ;
00883     if (subCode[n]==0 && allAnalVars2->getSize()>0) {
00884       coutE(InputArguments) << "RooAddPdf::getAnalyticalIntegral(" << GetName() << ") WARNING: component PDF " << pdf->GetName() 
00885                             << "   advertises inconsistent set of integrals (e.g. (X,Y) but not X or Y individually."
00886                             << "   Distributed analytical integration disabled. Please fix PDF" << endl ;
00887       allOK = kFALSE ;
00888     }
00889     delete allAnalVars2 ; 
00890     n++ ;
00891   }  
00892   if (!allOK) {
00893     delete[] subCode ;
00894     delete avIter ;
00895     return 0 ;
00896   }
00897 
00898   // Mare all analytically integrated observables as such
00899   analVars.add(allAnalVars) ;
00900 
00901   // Store set of variables analytically integrated
00902   RooArgSet* intSet = new RooArgSet(allAnalVars) ;
00903   Int_t masterCode = _codeReg.store(subCode,_pdfList.getSize(),intSet)+1 ;
00904 
00905   delete[] subCode ;
00906   delete avIter ;
00907 
00908   return masterCode ;
00909 }
00910 
00911 
00912 
00913 //_____________________________________________________________________________
00914 Double_t RooAddPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const 
00915 {
00916   // Return analytical integral defined by given scenario code
00917 
00918   // WVE needs adaptation to handle new rangeName feature
00919   if (code==0) {
00920     return getVal(normSet) ;
00921   }
00922 
00923   // Retrieve analytical integration subCodes and set of observabels integrated over
00924   RooArgSet* intSet ;
00925   const Int_t* subCode = _codeReg.retrieve(code-1,intSet) ;
00926   if (!subCode) {
00927     coutE(InputArguments) << "RooAddPdf::analyticalIntegral(" << GetName() << "): ERROR unrecognized integration code, " << code << endl ;
00928     assert(0) ;    
00929   }
00930 
00931   cxcoutD(Caching) << "RooAddPdf::aiWN(" << GetName() << ") calling getProjCache with nset = " << (normSet?*normSet:RooArgSet()) << endl ;
00932 
00933   if ((normSet==0 || normSet->getSize()==0) && _refCoefNorm.getSize()>0) {
00934 //     cout << "WVE integration of RooAddPdf without normalization, but have reference set, using ref set for normalization" << endl ;
00935     normSet = &_refCoefNorm ;
00936   }
00937 
00938   CacheElem* cache = getProjCache(normSet,intSet,0) ; // WVE rangename here?
00939   updateCoefficients(*cache,normSet) ;
00940 
00941   // Calculate the current value of this object  
00942   Double_t value(0) ;
00943 
00944   // Do running sum of coef/pdf pairs, calculate lastCoef.
00945   _pdfIter->Reset() ;
00946   _coefIter->Reset() ;
00947   RooAbsPdf* pdf ;
00948   Double_t snormVal ;
00949   Int_t i(0) ;
00950 
00951   //cout << "ROP::aIWN updateCoefCache with rangeName = " << (rangeName?rangeName:"<null>") << endl ;
00952   RooArgList* snormSet = (cache->_suppNormList.getSize()>0) ? &cache->_suppNormList : 0 ;
00953   while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
00954     if (_coefCache[i]) {
00955       snormVal = snormSet ? ((RooAbsReal*) cache->_suppNormList.at(i))->getVal() : 1.0 ;
00956       
00957       // WVE swap this?
00958       Double_t val = pdf->analyticalIntegralWN(subCode[i],normSet,rangeName) ;
00959       if (pdf->isSelectedComp()) {
00960         
00961         value += val*_coefCache[i]/snormVal ;
00962       }
00963     }    
00964     i++ ;
00965   }
00966 
00967   return value ;
00968 }
00969 
00970 
00971 
00972 //_____________________________________________________________________________
00973 Double_t RooAddPdf::expectedEvents(const RooArgSet* nset) const 
00974 {  
00975   // Return the number of expected events, which is either the sum of all coefficients
00976   // or the sum of the components extended terms, multiplied with the fraction that
00977   // is in the current range w.r.t the reference range
00978 
00979   Double_t expectedTotal(0.0);
00980 
00981   cxcoutD(Caching) << "RooAddPdf::expectedEvents(" << GetName() << ") calling getProjCache with nset = " << (nset?*nset:RooArgSet()) << endl ;
00982   CacheElem* cache = getProjCache(nset) ;
00983   updateCoefficients(*cache,nset) ;
00984 
00985   Int_t i(0) ;
00986   for (i=0 ; i<_pdfList.getSize() ; i++) {
00987 
00988     Double_t proj(1) ;
00989     RooAbsReal* r1 = ((RooAbsReal*)cache->_refRangeProjList.at(i)) ;
00990     RooAbsReal* r2 = ((RooAbsReal*)cache->_rangeProjList.at(i)) ;    
00991     if (r1 && r2) {
00992       proj = (r2->getVal()/r1->getVal()) ;      
00993     }
00994 
00995     Double_t ncomp =  _allExtendable ? ((RooAbsPdf*)_pdfList.at(i))->expectedEvents(nset) :  ((RooAbsReal*)_coefList.at(i))->getVal(nset) ;
00996     expectedTotal += proj*ncomp ;
00997   }
00998   return expectedTotal ;  
00999 }
01000 
01001 
01002 
01003 //_____________________________________________________________________________
01004 void RooAddPdf::selectNormalization(const RooArgSet* depSet, Bool_t force) 
01005 {
01006   // Interface function used by test statistics to freeze choice of observables
01007   // for interpretation of fraction coefficients
01008 
01009 
01010   if (!force && _refCoefNorm.getSize()!=0) {
01011     return ;
01012   }
01013 
01014   if (!depSet) {
01015     fixCoefNormalization(RooArgSet()) ;
01016     return ;
01017   }
01018 
01019   RooArgSet* myDepSet = getObservables(depSet) ;
01020   fixCoefNormalization(*myDepSet) ;
01021   delete myDepSet ;
01022 }
01023 
01024 
01025 
01026 //_____________________________________________________________________________
01027 void RooAddPdf::selectNormalizationRange(const char* rangeName, Bool_t force) 
01028 {
01029   // Interface function used by test statistics to freeze choice of range
01030   // for interpretation of fraction coefficients
01031 
01032   if (!force && _refCoefRangeName) {
01033     return ;
01034   }
01035 
01036   fixCoefRange(rangeName) ;
01037 }
01038 
01039 
01040 
01041 //_____________________________________________________________________________
01042 RooAbsGenContext* RooAddPdf::genContext(const RooArgSet &vars, const RooDataSet *prototype, 
01043                                         const RooArgSet* auxProto, Bool_t verbose) const 
01044 {
01045   // Return specialized context to efficiently generate toy events from RooAddPdfs
01046   // return RooAbsPdf::genContext(vars,prototype,auxProto,verbose) ; // WVE DEBUG
01047 
01048   return new RooAddGenContext(*this,vars,prototype,auxProto,verbose) ;
01049 }
01050 
01051 
01052 
01053 //_____________________________________________________________________________
01054 RooArgList RooAddPdf::CacheElem::containedArgs(Action) 
01055 {
01056   // List all RooAbsArg derived contents in this cache element
01057 
01058   RooArgList allNodes;
01059   allNodes.add(_projList) ;
01060   allNodes.add(_suppProjList) ;
01061   allNodes.add(_refRangeProjList) ;
01062   allNodes.add(_rangeProjList) ;
01063 
01064   return allNodes ;
01065 }
01066 
01067 
01068 
01069 //_____________________________________________________________________________
01070 std::list<Double_t>* RooAddPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const 
01071 {
01072   // Loop over components for plot sampling hints and merge them if there are multiple
01073 
01074   list<Double_t>* sumHint = 0 ;
01075 
01076   _pdfIter->Reset() ;
01077   RooAbsPdf* pdf ;
01078   // Loop over components pdf
01079   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
01080 
01081     list<Double_t>* pdfHint = pdf->plotSamplingHint(obs,xlo,xhi) ;
01082 
01083     // Process hint
01084     if (pdfHint) {
01085       if (!sumHint) {
01086 
01087         // If this is the first hint, then just save it
01088         sumHint = pdfHint ;
01089 
01090       } else {
01091 
01092         list<Double_t>* newSumHint = new list<Double_t>(sumHint->size()+pdfHint->size()) ;
01093 
01094         // Merge hints into temporary array
01095         merge(pdfHint->begin(),pdfHint->end(),sumHint->begin(),sumHint->end(),newSumHint->begin()) ;
01096 
01097         // Copy merged array without duplicates to new sumHintArrau
01098         delete sumHint ;
01099         sumHint = newSumHint ;
01100         
01101       }
01102     }
01103   }
01104 
01105   return sumHint ;
01106 }
01107 
01108 
01109 //_____________________________________________________________________________
01110 void RooAddPdf::printMetaArgs(ostream& os) const 
01111 {
01112   // Customized printing of arguments of a RooAddPdf to more intuitively reflect the contents of the
01113   // product operator construction
01114 
01115   _pdfIter->Reset() ;
01116   _coefIter->Reset() ;
01117 
01118   Bool_t first(kTRUE) ;
01119     
01120   RooAbsArg* coef, *pdf ;
01121   if (_coefList.getSize()!=0) { 
01122     while((coef=(RooAbsArg*)_coefIter->Next())) {
01123       if (!first) {
01124         os << " + " ;
01125       } else {
01126         first = kFALSE ;
01127       }
01128       pdf=(RooAbsArg*)_pdfIter->Next() ;
01129       os << coef->GetName() << " * " << pdf->GetName() ;
01130     }
01131     pdf = (RooAbsArg*) _pdfIter->Next() ;
01132     if (pdf) {
01133       os << " + [%] * " << pdf->GetName() ;
01134     }
01135   } else {
01136     
01137     while((pdf=(RooAbsArg*)_pdfIter->Next())) {
01138       if (!first) {
01139         os << " + " ;
01140       } else {
01141         first = kFALSE ;
01142       }
01143       os << pdf->GetName() ; 
01144     }  
01145   }
01146 
01147   os << " " ;    
01148 }

Generated on Tue Jul 5 15:06:00 2011 for ROOT_528-00b_version by  doxygen 1.5.1