RooAddModel.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooAddModel.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 // RooAddModel 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 // RooAddPdf relies on each component PDF to be normalized and will perform 
00037 // no normalization other than calculating the proper last coefficient c_n, if requested.
00038 // An (enforced) condition for this assuption is that each PDF_i is independent
00039 // of each coefficient_i.
00040 //
00041 // 
00042 
00043 #include "RooFit.h"
00044 #include "RooMsgService.h"
00045 
00046 #include "TIterator.h"
00047 #include "TIterator.h"
00048 #include "TList.h"
00049 #include "RooAddModel.h"
00050 #include "RooDataSet.h"
00051 #include "RooRealProxy.h"
00052 #include "RooPlot.h"
00053 #include "RooRealVar.h"
00054 #include "RooAddGenContext.h"
00055 #include "RooRealConstant.h"
00056 #include "RooNameReg.h"
00057 #include "RooMsgService.h"
00058 #include "RooRealIntegral.h"
00059 
00060 #include "Riostream.h"
00061 
00062 
00063 
00064 ClassImp(RooAddModel)
00065 ;
00066 
00067 
00068 //_____________________________________________________________________________
00069 RooAddModel::RooAddModel() :
00070   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00071   _refCoefRangeName(0),
00072   _codeReg(10),
00073   _snormList(0)
00074 {
00075   _pdfIter   = _pdfList.createIterator() ;
00076   _coefIter  = _coefList.createIterator() ;
00077 
00078   _coefCache = new Double_t[10] ;
00079   _coefErrCount = _errorCount ;
00080 }
00081 
00082 
00083 
00084 //_____________________________________________________________________________
00085 RooAddModel::RooAddModel(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, Bool_t ownPdfList) :
00086   RooResolutionModel(name,title,((RooResolutionModel*)inPdfList.at(0))->convVar()),
00087   _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
00088   _refCoefRangeName(0),
00089   _projectCoefs(kFALSE),
00090   _projCacheMgr(this,10),
00091   _intCacheMgr(this,10),
00092   _codeReg(10),
00093   _pdfList("!pdfs","List of PDFs",this),
00094   _coefList("!coefficients","List of coefficients",this),
00095   _haveLastCoef(kFALSE),
00096   _allExtendable(kFALSE)
00097 { 
00098   // Generic constructor from list of PDFs and list of coefficients.
00099   // Each pdf list element (i) is paired with coefficient list element (i).
00100   // The number of coefficients must be either equal to the number of PDFs,
00101   // in which case extended MLL fitting is enabled, or be one less.
00102   //
00103   // All PDFs must inherit from RooAbsPdf. All coefficients must inherit from RooAbsReal
00104 
00105   if (inPdfList.getSize()>inCoefList.getSize()+1) {
00106     coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() 
00107                           << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
00108     assert(0) ;
00109   }
00110 
00111   _pdfIter  = _pdfList.createIterator() ;
00112   _coefIter = _coefList.createIterator() ;
00113  
00114   // Constructor with N PDFs and N or N-1 coefs
00115   TIterator* pdfIter = inPdfList.createIterator() ;
00116   TIterator* coefIter = inCoefList.createIterator() ;
00117   RooAbsPdf* pdf ;
00118   RooAbsReal* coef ;
00119 
00120   while((coef = (RooAbsPdf*)coefIter->Next())) {
00121     pdf = (RooAbsPdf*) pdfIter->Next() ;
00122     if (!pdf) {
00123       coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() 
00124                             << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
00125       assert(0) ;
00126     }
00127     if (!dynamic_cast<RooAbsReal*>(coef)) {
00128       coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal, ignored" << endl ;
00129       continue ;
00130     }
00131     if (!dynamic_cast<RooAbsReal*>(pdf)) {
00132       coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") pdf " << pdf->GetName() << " is not of type RooAbsPdf, ignored" << endl ;
00133       continue ;
00134     }
00135     _pdfList.add(*pdf) ;
00136     _coefList.add(*coef) ;    
00137   }
00138 
00139   pdf = (RooAbsPdf*) pdfIter->Next() ;
00140   if (pdf) {
00141     if (!dynamic_cast<RooAbsReal*>(pdf)) {
00142       coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") last pdf " << coef->GetName() << " is not of type RooAbsPdf, fatal error" << endl ;
00143       assert(0) ;
00144     }
00145     _pdfList.add(*pdf) ;  
00146   } else {
00147     _haveLastCoef=kTRUE ;
00148   }
00149 
00150   delete pdfIter ;
00151   delete coefIter  ;
00152 
00153   _coefCache = new Double_t[_pdfList.getSize()] ;
00154   _coefErrCount = _errorCount ;
00155 
00156   if (ownPdfList) {
00157     _ownedComps.addOwned(_pdfList) ;
00158   }
00159 
00160 }
00161 
00162 
00163 
00164 //_____________________________________________________________________________
00165 RooAddModel::RooAddModel(const RooAddModel& other, const char* name) :
00166   RooResolutionModel(other,name),
00167   _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
00168   _refCoefRangeName((TNamed*)other._refCoefRangeName),
00169   _projectCoefs(other._projectCoefs),
00170   _projCacheMgr(other._projCacheMgr,this),
00171   _intCacheMgr(other._intCacheMgr,this),
00172   _codeReg(other._codeReg),
00173   _pdfList("!pdfs",this,other._pdfList),
00174   _coefList("!coefficients",this,other._coefList),
00175   _haveLastCoef(other._haveLastCoef),
00176   _allExtendable(other._allExtendable)
00177 {
00178   // Copy constructor
00179 
00180   _pdfIter  = _pdfList.createIterator() ;
00181   _coefIter = _coefList.createIterator() ;
00182   _coefCache = new Double_t[_pdfList.getSize()] ;
00183   _coefErrCount = _errorCount ;
00184 }
00185 
00186 
00187 
00188 //_____________________________________________________________________________
00189 RooAddModel::~RooAddModel()
00190 {
00191   // Destructor
00192 
00193   delete _pdfIter ;
00194   delete _coefIter ;
00195 
00196   if (_coefCache) delete[] _coefCache ;
00197 }
00198 
00199 
00200 
00201 //_____________________________________________________________________________
00202 void RooAddModel::fixCoefNormalization(const RooArgSet& refCoefNorm) 
00203 {
00204   // By default the interpretation of the fraction coefficients is
00205   // performed in the contextual choice of observables. This makes the
00206   // shape of the p.d.f explicitly dependent on the choice of
00207   // observables. This method instructs RooAddPdf to freeze the
00208   // interpretation of the coefficients to be done in the given set of
00209   // observables. If frozen, fractions are automatically transformed
00210   // from the reference normalization set to the contextual normalization
00211   // set by ratios of integrals
00212 
00213   if (refCoefNorm.getSize()==0) {
00214     _projectCoefs = kFALSE ;
00215     return ;
00216   }
00217   _projectCoefs = kTRUE ;  
00218 
00219   _refCoefNorm.removeAll() ;
00220   _refCoefNorm.add(refCoefNorm) ;
00221 
00222   _projCacheMgr.reset() ;
00223 }
00224 
00225 
00226 
00227 //_____________________________________________________________________________
00228 void RooAddModel::fixCoefRange(const char* rangeName)
00229 {
00230   // By default the interpretation of the fraction coefficients is
00231   // performed in the default range. This make the shape of a RooAddPdf
00232   // explicitly dependent on the range of the observables. To allow
00233   // a range independent definition of the fraction this function
00234   // instructs RooAddPdf to freeze its interpretation in the given
00235   // named range. If the current normalization range is different
00236   // from the reference range, the appropriate fraction coefficients
00237   // are automically calculation from the reference fractions using
00238   // ratios if integrals
00239 
00240   _refCoefRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
00241   if (_refCoefRangeName) _projectCoefs = kTRUE ;
00242 }
00243 
00244 
00245 
00246 //_____________________________________________________________________________
00247 RooResolutionModel* RooAddModel::convolution(RooFormulaVar* inBasis, RooAbsArg* owner) const
00248 {
00249   // Instantiate a clone of this resolution model representing a convolution with given
00250   // basis function. The owners object name is incorporated in the clones name
00251   // to avoid multiple convolution objects with the same name in complex PDF structures.
00252   //
00253   // RooAddModel will clone all the component models to create a composite convolution object
00254 
00255   // Check that primary variable of basis functions is our convolution variable  
00256   if (inBasis->findServer(0) != x.absArg()) {
00257     coutE(InputArguments) << "RooAddModel::convolution(" << GetName() 
00258                           << ") convolution parameter of basis function and PDF don't match" << endl ;
00259     ccoutE(InputArguments) << "basis->findServer(0) = " << inBasis->findServer(0) << " " << inBasis->findServer(0)->GetName() << endl ;
00260     ccoutE(InputArguments) << "x.absArg()           = " << x.absArg() << " " << x.absArg()->GetName() << endl ;
00261     inBasis->Print("v") ;
00262     return 0 ;
00263   }
00264 
00265   TString newName(GetName()) ;
00266   newName.Append("_conv_") ;
00267   newName.Append(inBasis->GetName()) ;
00268   newName.Append("_[") ;
00269   newName.Append(owner->GetName()) ;
00270   newName.Append("]") ;
00271 
00272   TString newTitle(GetTitle()) ;
00273   newTitle.Append(" convoluted with basis function ") ;
00274   newTitle.Append(inBasis->GetName()) ;
00275 
00276   _pdfIter->Reset() ;
00277   RooResolutionModel* model ;
00278   RooArgList modelList ;
00279   while((model = (RooResolutionModel*)_pdfIter->Next())) {       
00280     // Create component convolution
00281     RooResolutionModel* conv = model->convolution(inBasis,owner) ;    
00282     modelList.add(*conv) ;
00283   }
00284 
00285   _coefIter->Reset() ;
00286   RooAbsReal* coef ;
00287   RooArgList theCoefList ;  
00288   while((coef = (RooAbsReal*)_coefIter->Next())) {
00289     theCoefList.add(*coef) ;
00290   }
00291     
00292   RooAddModel* convSum = new RooAddModel(newName,newTitle,modelList,theCoefList,kTRUE) ;
00293   convSum->changeBasis(inBasis) ;
00294   return convSum ;
00295 }
00296 
00297 
00298 
00299 //_____________________________________________________________________________
00300 Int_t RooAddModel::basisCode(const char* name) const 
00301 {
00302   // Return code for basis function representing by 'name' string.
00303   // The basis code of the first component model will be returned,
00304   // if the basis is supported by all components. Otherwise 0
00305   // is returned
00306 
00307   TIterator* mIter = _pdfList.createIterator() ;
00308   RooResolutionModel* model ;
00309   Bool_t first(kTRUE), code(0) ;
00310     while((model = (RooResolutionModel*)mIter->Next())) {
00311       Int_t subCode = model->basisCode(name) ;
00312       if (first) {
00313         code = subCode ;
00314         first = kFALSE ;
00315       } else if (subCode==0) {
00316         code = 0 ;
00317       }
00318   }
00319   delete mIter ;
00320 
00321   return code ;
00322 }
00323 
00324 
00325 
00326 //_____________________________________________________________________________
00327 RooAddModel::CacheElem* RooAddModel::getProjCache(const RooArgSet* nset, const RooArgSet* iset, const char* rangeName) const
00328 {
00329   // Retrieve cache element with for calculation of p.d.f value with normalization set nset and integrated over iset
00330   // in range 'rangeName'. If cache element does not exist, create and fill it on the fly. The cache contains
00331   // suplemental normalization terms (in case not all added p.d.f.s have the same observables), projection
00332   // integrals to calculated transformed fraction coefficients when a frozen reference frame is provided
00333   // and projection integrals for similar transformations when a frozen reference range is provided.
00334 
00335   // Check if cache already exists 
00336   CacheElem* cache = (CacheElem*) _projCacheMgr.getObj(nset,iset,0,RooNameReg::ptr(rangeName)) ;
00337   if (cache) {
00338     return cache ;
00339   }
00340 
00341   //Create new cache 
00342   cache = new CacheElem ;
00343 
00344   // *** PART 1 : Create supplemental normalization list ***
00345 
00346   // Retrieve the combined set of dependents of this PDF ;
00347   RooArgSet *fullDepList = getObservables(nset) ;
00348   if (iset) {
00349     fullDepList->remove(*iset,kTRUE,kTRUE) ;
00350   }    
00351 
00352   // Fill with dummy unit RRVs for now
00353   _pdfIter->Reset() ;
00354   _coefIter->Reset() ;
00355   RooAbsPdf* pdf ;
00356   RooAbsReal* coef ;
00357   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {    
00358     coef=(RooAbsPdf*)_coefIter->Next() ;
00359 
00360     // Start with full list of dependents
00361     RooArgSet supNSet(*fullDepList) ;
00362 
00363     // Remove PDF dependents
00364     RooArgSet* pdfDeps = pdf->getObservables(nset) ;
00365     if (pdfDeps) {
00366       supNSet.remove(*pdfDeps,kTRUE,kTRUE) ;
00367       delete pdfDeps ; 
00368     }
00369 
00370     // Remove coef dependents
00371     RooArgSet* coefDeps = coef ? coef->getObservables(nset) : 0 ;
00372     if (coefDeps) {
00373       supNSet.remove(*coefDeps,kTRUE,kTRUE) ;
00374       delete coefDeps ;
00375     }
00376     
00377     RooAbsReal* snorm ;
00378     TString name(GetName()) ;
00379     name.Append("_") ;
00380     name.Append(pdf->GetName()) ;
00381     name.Append("_SupNorm") ;
00382     if (supNSet.getSize()>0) {
00383       snorm = new RooRealIntegral(name,"Supplemental normalization integral",RooRealConstant::value(1.0),supNSet) ;
00384     } else {
00385       snorm = new RooRealVar(name,"Unit Supplemental normalization integral",1.0) ;
00386     }
00387     cache->_suppNormList.addOwned(*snorm) ;
00388   }
00389 
00390   delete fullDepList ;
00391     
00392   if (_verboseEval>1) {
00393     cxcoutD(Caching) << "RooAddModel::syncSuppNormList(" << GetName() << ") synching supplemental normalization list for norm" << (nset?*nset:RooArgSet()) << endl ;
00394     if (dologD(Caching)) {
00395       cache->_suppNormList.Print("v") ;
00396     }
00397   }
00398 
00399 
00400   // *** PART 2 : Create projection coefficients ***
00401 
00402   // If no projections required stop here
00403   if (!_projectCoefs || _basis!=0 ) {
00404     _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
00405     return cache ;
00406   }
00407 
00408 
00409   // Reduce iset/nset to actual dependents of this PDF
00410   RooArgSet* nset2 = nset ? getObservables(nset) : new RooArgSet() ;
00411 
00412   // Check if requested transformation is not identity 
00413   if (!nset2->equals(_refCoefNorm) || _refCoefRangeName !=0 || rangeName !=0 ) {
00414     
00415     coutI(Caching)  << "RooAddModel::syncCoefProjList(" << GetName() << ") creating coefficient projection integrals" << endl ;
00416     ccoutI(Caching) << "  from current normalization: "  ; nset2->Print("1") ;
00417     ccoutI(Caching) << "          with current range: " << (rangeName?rangeName:"<none>") << endl ;
00418     ccoutI(Caching) << "  to reference normalization: "  ; _refCoefNorm.Print("1") ; 
00419     ccoutI(Caching) << "        with reference range: " << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"<none>") << endl ;
00420     
00421     // Recalculate projection integrals of PDFs 
00422     _pdfIter->Reset() ;
00423     RooAbsPdf* thePdf ;
00424 
00425     while((thePdf=(RooAbsPdf*)_pdfIter->Next())) {
00426 
00427       // Calculate projection integral
00428       RooAbsReal* pdfProj ;
00429       if (!nset2->equals(_refCoefNorm)) {
00430         pdfProj = thePdf->createIntegral(*nset2,_refCoefNorm) ;
00431         pdfProj->setOperMode(operMode()) ;
00432       } else {
00433         TString name(GetName()) ;
00434         name.Append("_") ;
00435         name.Append(thePdf->GetName()) ;
00436         name.Append("_ProjectNorm") ;
00437         pdfProj = new RooRealVar(name,"Unit Projection normalization integral",1.0) ;
00438       }
00439 
00440       cache->_projList.addOwned(*pdfProj) ;
00441 
00442       // Calculation optional supplemental normalization term
00443       RooArgSet supNormSet(_refCoefNorm) ;
00444       RooArgSet* deps = thePdf->getParameters(RooArgSet()) ;
00445       supNormSet.remove(*deps,kTRUE,kTRUE) ;
00446       delete deps ;
00447 
00448       RooAbsReal* snorm ;
00449       TString name(GetName()) ;
00450       name.Append("_") ;
00451       name.Append(thePdf->GetName()) ;
00452       name.Append("_ProjSupNorm") ;
00453       if (supNormSet.getSize()>0) {
00454         snorm = new RooRealIntegral(name,"Projection Supplemental normalization integral",
00455                                     RooRealConstant::value(1.0),supNormSet) ;
00456       } else {
00457         snorm = new RooRealVar(name,"Unit Projection Supplemental normalization integral",1.0) ;
00458       }
00459       cache->_suppProjList.addOwned(*snorm) ;
00460 
00461       // Calculate reference range adjusted projection integral
00462       RooAbsReal* rangeProj1 ;
00463       if (_refCoefRangeName && _refCoefNorm.getSize()>0) {
00464         rangeProj1 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,RooNameReg::str(_refCoefRangeName)) ;
00465         rangeProj1->setOperMode(operMode()) ;
00466       } else {
00467         TString theName(GetName()) ;
00468         theName.Append("_") ;
00469         theName.Append(thePdf->GetName()) ;
00470         theName.Append("_RangeNorm1") ;
00471         rangeProj1 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
00472       }
00473       cache->_refRangeProjList.addOwned(*rangeProj1) ;
00474       
00475 
00476       // Calculate range adjusted projection integral
00477       RooAbsReal* rangeProj2 ;
00478       if (rangeName && _refCoefNorm.getSize()>0) {
00479         rangeProj2 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,rangeName) ;
00480         rangeProj2->setOperMode(operMode()) ;
00481       } else {
00482         TString theName(GetName()) ;
00483         theName.Append("_") ;
00484         theName.Append(thePdf->GetName()) ;
00485         theName.Append("_RangeNorm2") ;
00486         rangeProj2 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
00487       }
00488       cache->_rangeProjList.addOwned(*rangeProj2) ;
00489 
00490     }               
00491 
00492   }
00493 
00494   delete nset2 ;
00495 
00496   _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
00497 
00498   return cache ;
00499 }
00500 
00501 
00502 
00503 //_____________________________________________________________________________
00504 void RooAddModel::updateCoefficients(CacheElem& cache, const RooArgSet* nset) const 
00505 {
00506   // Update the coefficient values in the given cache element: calculate new remainder
00507   // fraction, normalize fractions obtained from extended ML terms to unity and
00508   // multiply these the various range and dimensional corrections needed in the
00509   // current use context
00510 
00511   // cxcoutD(ChangeTracking) << "RooAddModel::updateCoefficients(" << GetName() << ") update coefficients" << endl ;
00512   
00513   Int_t i ;
00514 
00515   // Straight coefficients
00516   if (_allExtendable) {
00517     
00518     // coef[i] = expectedEvents[i] / SUM(expectedEvents)
00519     Double_t coefSum(0) ;
00520     for (i=0 ; i<_pdfList.getSize() ; i++) {
00521       _coefCache[i] = ((RooAbsPdf*)_pdfList.at(i))->expectedEvents(_refCoefNorm.getSize()>0?&_refCoefNorm:nset) ;
00522       coefSum += _coefCache[i] ;
00523     }
00524     if (coefSum==0.) {
00525       coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName() << ") WARNING: total number of expected events is 0" << endl ;
00526     } else {
00527       for (i=0 ; i<_pdfList.getSize() ; i++) {
00528         _coefCache[i] /= coefSum ;
00529       }                             
00530     }
00531     
00532   } else {
00533     if (_haveLastCoef) {
00534       
00535       // coef[i] = coef[i] / SUM(coef)
00536       Double_t coefSum(0) ;
00537       for (i=0 ; i<_coefList.getSize() ; i++) {
00538         _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
00539         coefSum += _coefCache[i] ;
00540       }         
00541       for (i=0 ; i<_coefList.getSize() ; i++) {
00542         _coefCache[i] /= coefSum ;
00543       }                 
00544     } else {
00545       
00546       // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
00547       Double_t lastCoef(1) ;
00548       for (i=0 ; i<_coefList.getSize() ; i++) {
00549         _coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
00550         cxcoutD(Caching) << "SYNC: orig coef[" << i << "] = " << _coefCache[i] << endl ;
00551         lastCoef -= _coefCache[i] ;
00552       }                 
00553       _coefCache[_coefList.getSize()] = lastCoef ;
00554       cxcoutD(Caching) << "SYNC: orig coef[" << _coefList.getSize() << "] = " << _coefCache[_coefList.getSize()] << endl ;
00555       
00556       
00557       // Warn about coefficient degeneration
00558       if ((lastCoef<-1e-05 || (lastCoef-1)>1e-5) && _coefErrCount-->0) {
00559         coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName() 
00560                     << " WARNING: sum of PDF coefficients not in range [0-1], value=" 
00561                     << 1-lastCoef << endl ;
00562         if (_coefErrCount==0) {
00563           coutW(Eval) << " (no more will be printed)" << endl  ;
00564         }
00565       } 
00566     }
00567   }
00568 
00569   
00570 
00571   // Stop here if not projection is required or needed
00572   if ((!_projectCoefs) || cache._projList.getSize()==0) {
00573     //     cout << "SYNC no projection required rangeName = " << (rangeName?rangeName:"<none>") << endl ;
00574     return ;
00575   }
00576 
00577   // Adjust coefficients for given projection
00578   Double_t coefSum(0) ;
00579   for (i=0 ; i<_pdfList.getSize() ; i++) {
00580     RooAbsPdf::globalSelectComp(kTRUE) ;    
00581 
00582     RooAbsReal* pp = ((RooAbsReal*)cache._projList.at(i)) ; 
00583     RooAbsReal* sn = ((RooAbsReal*)cache._suppProjList.at(i)) ; 
00584     RooAbsReal* r1 = ((RooAbsReal*)cache._refRangeProjList.at(i)) ;
00585     RooAbsReal* r2 = ((RooAbsReal*)cache._rangeProjList.at(i)) ;
00586     
00587     if (dologD(Eval)) {
00588       cxcoutD(Eval) << "pp = " << pp->GetName() << endl 
00589                     << "sn = " << sn->GetName() << endl 
00590                     << "r1 = " << r1->GetName() << endl 
00591                     << "r2 = " << r2->GetName() << endl ;
00592       r1->printStream(ccoutD(Eval),kName|kArgs|kValue,kSingleLine) ;
00593       r1->printCompactTree(ccoutD(Eval)) ;
00594     }
00595 
00596     Double_t proj = pp->getVal()/sn->getVal()*(r2->getVal()/r1->getVal()) ;  
00597     
00598     RooAbsPdf::globalSelectComp(kFALSE) ;
00599 
00600     _coefCache[i] *= proj ;
00601     coefSum += _coefCache[i] ;
00602   }
00603   for (i=0 ; i<_pdfList.getSize() ; i++) {
00604     _coefCache[i] /= coefSum ;
00605 //     cout << "POST-SYNC coef[" << i << "] = " << _coefCache[i] << endl ;
00606   }
00607    
00608 
00609   
00610 }
00611 
00612 
00613 
00614 //_____________________________________________________________________________
00615 Double_t RooAddModel::evaluate() const 
00616 {
00617   // Calculate the current value
00618   const RooArgSet* nset = _normSet ; 
00619   CacheElem* cache = getProjCache(nset) ;
00620 
00621   updateCoefficients(*cache,nset) ;
00622 
00623   
00624   // Do running sum of coef/pdf pairs, calculate lastCoef.
00625   _pdfIter->Reset() ;
00626   _coefIter->Reset() ;
00627   RooAbsPdf* pdf ;
00628 
00629   Double_t snormVal ;
00630   Double_t value(0) ;
00631   Int_t i(0) ;
00632   while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
00633     if (_coefCache[i]!=0.) {
00634       snormVal = nset ? ((RooAbsReal*)cache->_suppNormList.at(i))->getVal() : 1.0 ;
00635       Double_t pdfVal = pdf->getVal(nset) ;
00636       // Double_t pdfNorm = pdf->getNorm(nset) ;
00637       if (pdf->isSelectedComp()) {
00638         value += pdfVal*_coefCache[i]/snormVal ;
00639         cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ")  value += [" 
00640                         << pdf->GetName() << "] " << pdfVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
00641       }
00642     }
00643     i++ ;
00644   }
00645 
00646   return value ;
00647 }
00648 
00649 
00650 
00651 //_____________________________________________________________________________
00652 void RooAddModel::resetErrorCounters(Int_t resetValue)
00653 {
00654   // Reset error counter to given value, limiting the number
00655   // of future error messages for this pdf to 'resetValue'
00656   RooAbsPdf::resetErrorCounters(resetValue) ;
00657   _coefErrCount = resetValue ;
00658 }
00659 
00660 
00661 
00662 //_____________________________________________________________________________
00663 Bool_t RooAddModel::checkObservables(const RooArgSet* nset) const 
00664 {
00665   // Check if PDF is valid for given normalization set.
00666   // Coeffient and PDF must be non-overlapping, but pdf-coefficient 
00667   // pairs may overlap each other
00668 
00669   Bool_t ret(kFALSE) ;
00670 
00671   _pdfIter->Reset() ;
00672   _coefIter->Reset() ;
00673   RooAbsReal* coef ;
00674   RooAbsReal* pdf ;
00675   while((coef=(RooAbsReal*)_coefIter->Next())) {
00676     pdf = (RooAbsReal*)_pdfIter->Next() ;
00677     if (pdf->observableOverlaps(nset,*coef)) {
00678       coutE(InputArguments) << "RooAddModel::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName() 
00679                             << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
00680       ret = kTRUE ;
00681     }
00682   }
00683   
00684   return ret ;
00685 }
00686 
00687 
00688 
00689 //_____________________________________________________________________________
00690 Int_t RooAddModel::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
00691                                          const RooArgSet* normSet, const char* rangeName) const 
00692 {
00693 
00694   if (_forceNumInt) return 0 ;
00695 
00696   // Declare that we can analytically integrate all requested observables
00697   analVars.add(allVars) ;
00698 
00699   // Retrieve (or create) the required component integral list
00700   Int_t code ;
00701   RooArgList *cilist ;
00702   getCompIntList(normSet,&allVars,cilist,code,rangeName) ;
00703   
00704   return code+1 ;
00705   
00706 }
00707 
00708 
00709 
00710 //_____________________________________________________________________________
00711 void RooAddModel::getCompIntList(const RooArgSet* nset, const RooArgSet* iset, pRooArgList& compIntList, Int_t& code, const char* isetRangeName) const 
00712 {
00713   // Check if this configuration was created before
00714   Int_t sterileIdx(-1) ;
00715 
00716   IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObj(nset,iset,&sterileIdx,RooNameReg::ptr(isetRangeName)) ;
00717   if (cache) {
00718     code = _intCacheMgr.lastIndex() ;
00719     compIntList = &cache->_intList ;
00720     
00721     return ;
00722   }
00723 
00724   // Create containers for partial integral components to be generated
00725   cache = new IntCacheElem ;
00726 
00727   // Fill Cache
00728   _pdfIter->Reset() ;
00729   RooResolutionModel* model ;
00730   while ((model=(RooResolutionModel*)_pdfIter->Next())) {
00731     RooAbsReal* intPdf = model->createIntegral(*iset,nset,0,isetRangeName) ;
00732     cache->_intList.addOwned(*intPdf) ;
00733   }
00734 
00735   // Store the partial integral list and return the assigned code ;
00736   code = _intCacheMgr.setObj(nset,iset,(RooAbsCacheElement*)cache,RooNameReg::ptr(isetRangeName)) ;
00737 
00738   // Fill references to be returned
00739   compIntList = &cache->_intList ;
00740 }
00741 
00742 
00743 
00744 //_____________________________________________________________________________
00745 Double_t RooAddModel::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const 
00746 {
00747   // Return analytical integral defined by given scenario code
00748 
00749   // No integration scenario
00750   if (code==0) {
00751     return getVal(normSet) ;
00752   }
00753 
00754   // Partial integration scenarios
00755   IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObjByIndex(code-1) ;
00756   
00757   RooArgList* compIntList ;
00758 
00759   // If cache has been sterilized, revive this slot
00760   if (cache==0) {
00761     RooArgSet* vars = getParameters(RooArgSet()) ;
00762     RooArgSet* nset = _intCacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
00763     RooArgSet* iset = _intCacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
00764 
00765     Int_t code2(-1) ;
00766     getCompIntList(nset,iset,compIntList,code2,rangeName) ;
00767 
00768     delete vars ;
00769     delete nset ;
00770     delete iset ;
00771   } else {
00772 
00773     compIntList = &cache->_intList ;
00774 
00775   }
00776 
00777   // Calculate the current value
00778   const RooArgSet* nset = _normSet ; 
00779   CacheElem* pcache = getProjCache(nset) ;
00780 
00781   updateCoefficients(*pcache,nset) ;
00782   
00783   // Do running sum of coef/pdf pairs, calculate lastCoef.
00784   TIterator* compIntIter = compIntList->createIterator() ;
00785   _coefIter->Reset() ;
00786   RooAbsReal* pdfInt ;
00787 
00788   Double_t snormVal ;
00789   Double_t value(0) ;
00790   Int_t i(0) ;
00791   while((pdfInt = (RooAbsReal*)compIntIter->Next())) {
00792     if (_coefCache[i]!=0.) {
00793       snormVal = nset ? ((RooAbsReal*)pcache->_suppNormList.at(i))->getVal() : 1.0 ;
00794       Double_t intVal = pdfInt->getVal(nset) ;
00795       value += intVal*_coefCache[i]/snormVal ;
00796       cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ")  value += [" 
00797                       << pdfInt->GetName() << "] " << intVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
00798     }
00799     i++ ;
00800   }
00801 
00802   delete compIntIter ;
00803   
00804   return value ;
00805   
00806 }
00807 
00808 
00809 
00810 //_____________________________________________________________________________
00811 Double_t RooAddModel::expectedEvents(const RooArgSet* nset) const 
00812 {  
00813   // Return the number of expected events, which is either the sum of all coefficients
00814   // or the sum of the components extended terms
00815 
00816   Double_t expectedTotal(0.0);
00817   RooAbsPdf* pdf ;
00818     
00819   if (_allExtendable) {
00820     
00821     // Sum of the extended terms
00822     _pdfIter->Reset() ;
00823     while((pdf = (RooAbsPdf*)_pdfIter->Next())) {      
00824       expectedTotal += pdf->expectedEvents(nset) ;
00825     }   
00826     
00827   } else {
00828     
00829     // Sum the coefficients
00830     _coefIter->Reset() ;
00831     RooAbsReal* coef ;
00832     while((coef=(RooAbsReal*)_coefIter->Next())) {
00833       expectedTotal += coef->getVal() ;
00834     }   
00835   }
00836 
00837   return expectedTotal;
00838 }
00839 
00840 
00841 
00842 //_____________________________________________________________________________
00843 void RooAddModel::selectNormalization(const RooArgSet* depSet, Bool_t force) 
00844 {
00845   // Interface function used by test statistics to freeze choice of observables
00846   // for interpretation of fraction coefficients
00847 
00848   if (!force && _refCoefNorm.getSize()!=0) {
00849     return ;
00850   }
00851 
00852   if (!depSet) {
00853     fixCoefNormalization(RooArgSet()) ;
00854     return ;
00855   }
00856 
00857   RooArgSet* myDepSet = getObservables(depSet) ;
00858   fixCoefNormalization(*myDepSet) ;
00859   delete myDepSet ;
00860 }
00861 
00862 
00863 
00864 //_____________________________________________________________________________
00865 void RooAddModel::selectNormalizationRange(const char* rangeName, Bool_t force) 
00866 {
00867   // Interface function used by test statistics to freeze choice of range
00868   // for interpretation of fraction coefficients
00869 
00870   if (!force && _refCoefRangeName) {
00871     return ;
00872   }
00873 
00874   fixCoefRange(rangeName) ;
00875 }
00876 
00877 
00878 
00879 //_____________________________________________________________________________
00880 RooAbsGenContext* RooAddModel::genContext(const RooArgSet &vars, const RooDataSet *prototype, 
00881                                         const RooArgSet* auxProto, Bool_t verbose) const 
00882 {
00883   // Return specialized context to efficiently generate toy events from RooAddPdfs
00884 
00885   return new RooAddGenContext(*this,vars,prototype,auxProto,verbose) ;
00886 }
00887 
00888 
00889 
00890 //_____________________________________________________________________________
00891 Bool_t RooAddModel::isDirectGenSafe(const RooAbsArg& arg) const 
00892 {
00893   // Direct generation is safe if all components say so
00894   _pdfIter->Reset() ;
00895   RooAbsPdf* pdf ;
00896   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
00897     if (!pdf->isDirectGenSafe(arg)) {
00898       return kFALSE ;
00899     }
00900   }
00901   return kTRUE ;
00902 }
00903 
00904 
00905 
00906 //_____________________________________________________________________________
00907 Int_t RooAddModel::getGenerator(const RooArgSet& directVars, RooArgSet &/*generateVars*/, Bool_t /*staticInitOK*/) const
00908 {
00909   // Return pseud-code that indicates if all components can do internal generation (1) or not (0)
00910 
00911   _pdfIter->Reset() ;
00912   RooAbsPdf* pdf ;
00913   while((pdf=(RooAbsPdf*)_pdfIter->Next())) {
00914     RooArgSet tmp ;
00915     if (pdf->getGenerator(directVars,tmp)==0) {
00916       return 0 ;
00917     }
00918   }
00919   return 1 ;  
00920 }
00921 
00922 
00923 
00924 
00925 //_____________________________________________________________________________
00926 void RooAddModel::generateEvent(Int_t /*code*/)
00927 {
00928   // This function should never be called as RooAddModel implements a custom generator context
00929   assert(0) ;
00930 }
00931 
00932 
00933 
00934 
00935 //_____________________________________________________________________________
00936 RooArgList RooAddModel::CacheElem::containedArgs(Action) 
00937 {
00938   // List all RooAbsArg derived contents in this cache element
00939 
00940   RooArgList allNodes;
00941   allNodes.add(_projList) ;
00942   allNodes.add(_suppProjList) ;
00943   allNodes.add(_refRangeProjList) ;
00944   allNodes.add(_rangeProjList) ;
00945 
00946   return allNodes ;
00947 }
00948 
00949 
00950 
00951 //_____________________________________________________________________________
00952 RooArgList RooAddModel::IntCacheElem::containedArgs(Action) 
00953 {
00954   // List all RooAbsArg derived contents in this cache element
00955 
00956   RooArgList allNodes(_intList) ;
00957   return allNodes ;
00958 }
00959 
00960 
00961 //_____________________________________________________________________________
00962 void RooAddModel::printMetaArgs(ostream& os) const 
00963 {
00964   // Customized printing of arguments of a RooAddModel to more intuitively reflect the contents of the
00965   // product operator construction
00966 
00967   _pdfIter->Reset() ;
00968   _coefIter->Reset() ;
00969 
00970   Bool_t first(kTRUE) ;
00971     
00972   os << "(" ;
00973   RooAbsArg* coef, *pdf ;
00974   while((coef=(RooAbsArg*)_coefIter->Next())) {
00975     if (!first) {
00976       os << " + " ;
00977     } else {
00978       first = kFALSE ;
00979     }
00980     pdf=(RooAbsArg*)_pdfIter->Next() ;
00981     os << coef->GetName() << " * " << pdf->GetName() ;
00982   }
00983   pdf = (RooAbsArg*) _pdfIter->Next() ;
00984   if (pdf) {
00985     os << " + [%] * " << pdf->GetName() ;
00986   }
00987   os << ") " ;    
00988 }
00989 

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