RooProjectedPdf.cxx

Go to the documentation of this file.
00001  /***************************************************************************** 
00002   * Project: RooFit                                                           * 
00003   *                                                                           * 
00004   * Copyright (c) 2000-2005, Regents of the University of California          * 
00005   *                          and Stanford University. All rights reserved.    * 
00006   *                                                                           * 
00007   * Redistribution and use in source and binary forms,                        * 
00008   * with or without modification, are permitted according to the terms        * 
00009   * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             * 
00010   *****************************************************************************/ 
00011 
00012 //////////////////////////////////////////////////////////////////////////////
00013 //
00014 // BEGIN_HTML
00015 // Class RooProjectedPdf is a RooAbsPdf implementation that represent a projection 
00016 // of a given input p.d.f and the object returned by RooAbsPdf::createProjection.
00017 // <p>
00018 // The actual projection integral for it value and normalization are
00019 // calculated on the fly in getVal() once the normalization observables are known.
00020 // Class RooProjectedPdf can cache projected p.d.f.s for multiple normalization
00021 // observables simultaneously.
00022 // <p>
00023 // The createProjection() method of RooProjectedPdf is overloaded and will
00024 // return a new RooProjectedPdf that performs the projection of itself
00025 // and the requested additional projections in one integration step
00026 // The performance of <pre>f->createProjection(x)->createProjection(y)</pre>
00027 // is therefore identical to that of <pre>f->createProjection(RooArgSet(x,y))</pre>
00028 // END_HTML
00029 //
00030 
00031 #include "Riostream.h" 
00032 
00033 #include "RooFit.h"
00034 #include "RooProjectedPdf.h" 
00035 #include "RooMsgService.h"
00036 #include "RooAbsReal.h" 
00037 #include "RooRealVar.h"
00038 #include "RooNameReg.h"
00039 
00040 
00041  ClassImp(RooProjectedPdf) 
00042    ;
00043 
00044 
00045 //_____________________________________________________________________________
00046 RooProjectedPdf::RooProjectedPdf() : _curNormSet(0)
00047 {
00048   // Default constructor
00049 }
00050 
00051 
00052 
00053 //_____________________________________________________________________________
00054  RooProjectedPdf::RooProjectedPdf(const char *name, const char *title, RooAbsReal& _intpdf, const RooArgSet& intObs) :
00055    RooAbsPdf(name,title), 
00056    intpdf("!IntegratedPdf","intpdf",this,_intpdf,kFALSE,kFALSE),
00057    intobs("!IntegrationObservables","intobs",this,kFALSE,kFALSE),
00058    deps("!Dependents","deps",this,kTRUE,kTRUE),
00059    _cacheMgr(this,10)
00060  { 
00061    // Construct projection of input pdf '_intpdf' over observables 'intObs'
00062 
00063    intobs.add(intObs) ;
00064 
00065    // Add all other dependens of projected p.d.f. directly
00066    RooArgSet* tmpdeps = _intpdf.getParameters(intObs) ;
00067    deps.add(*tmpdeps) ;
00068    delete tmpdeps ;
00069  } 
00070 
00071 
00072 
00073 //_____________________________________________________________________________
00074  RooProjectedPdf::RooProjectedPdf(const RooProjectedPdf& other, const char* name) :  
00075    RooAbsPdf(other,name), 
00076    intpdf("!IntegratedPdf",this,other.intpdf),
00077    intobs("!IntegrationObservable",this,other.intobs),
00078    deps("!Dependents",this,other.deps),
00079    _cacheMgr(other._cacheMgr,this)
00080 { 
00081    // Copy constructor
00082  } 
00083 
00084 
00085 
00086 //_____________________________________________________________________________
00087 Double_t RooProjectedPdf::getVal(const RooArgSet* set) const 
00088 {
00089   // Special version of getVal() overrides RooAbsReal::getVal() to save value of current normalization set
00090 
00091   _curNormSet = (RooArgSet*)set ;
00092 
00093   return RooAbsPdf::getVal(set) ;
00094 }
00095 
00096 
00097 
00098 //_____________________________________________________________________________
00099 Double_t RooProjectedPdf::evaluate() const 
00100 {
00101   // Evaluate projected p.d.f
00102 
00103   // Calculate current unnormalized value of object
00104   int code ;
00105   const RooAbsReal* proj = getProjection(&intobs,_curNormSet,0,code) ;
00106   
00107   return proj->getVal() ;
00108 }
00109 
00110 
00111 
00112 //_____________________________________________________________________________
00113 const RooAbsReal* RooProjectedPdf::getProjection(const RooArgSet* iset, const RooArgSet* nset, const char* rangeName, int& code) const
00114 {
00115   // Retrieve object representing projection integral of input p.d.f
00116   // over observables iset, while normalizing over observables
00117   // nset. The code argument returned by reference is the unique code
00118   // defining this particular projection configuration
00119 
00120 
00121   // Check if this configuration was created before
00122   Int_t sterileIdx(-1) ;
00123   CacheElem* cache = (CacheElem*) _cacheMgr.getObj(iset,nset,&sterileIdx,RooNameReg::ptr(rangeName)) ;
00124   if (cache) {
00125     code = _cacheMgr.lastIndex() ;
00126     return static_cast<const RooAbsReal*>(cache->_projection);
00127   }
00128 
00129   RooArgSet* nset2 =  intpdf.arg().getObservables(*nset) ;
00130 
00131   if (iset) {
00132     nset2->add(*iset) ;
00133   }
00134   RooAbsReal* proj = intpdf.arg().createIntegral(iset?*iset:RooArgSet(),nset2,0,rangeName) ;
00135   delete nset2 ;
00136 
00137   cache = new CacheElem ;
00138   cache->_projection = proj ;
00139 
00140   code = _cacheMgr.setObj(iset,nset,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
00141 
00142   coutI(Integration) << "RooProjectedPdf::getProjection(" << GetName() << ") creating new projection " << proj->GetName() << " with code " << code << endl ;
00143 
00144   return proj ;
00145 }
00146 
00147 
00148 
00149 //_____________________________________________________________________________
00150 RooAbsPdf* RooProjectedPdf::createProjection(const RooArgSet& iset) 
00151 {
00152   // Special version of RooAbsReal::createProjection that deals with
00153   // projections of projections. Instead of integrating twice, a new
00154   // RooProjectedPdf is returned that is configured to perform the
00155   // complete integration in one step
00156 
00157   RooArgSet combiset(iset) ;
00158   combiset.add(intobs) ;
00159   return static_cast<RooAbsPdf&>( const_cast<RooAbsReal&>(intpdf.arg()) ).createProjection(combiset) ;
00160 }
00161 
00162 
00163 
00164 //_____________________________________________________________________________
00165 Bool_t RooProjectedPdf::forceAnalyticalInt(const RooAbsArg& /*dep*/) const 
00166 {
00167   // Force RooRealIntegral to relegate integration of all observables to internal logic
00168 
00169   return kTRUE ;
00170 }
00171 
00172 
00173 
00174 //_____________________________________________________________________________
00175 Int_t RooProjectedPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName) const 
00176 { 
00177   // Mark all requested variables as internally integrated
00178 
00179   analVars.add(allVars) ;
00180   
00181   // Create the appropriate integral
00182   int code ;
00183   RooArgSet allVars2(allVars) ;
00184   allVars2.add(intobs) ;
00185   getProjection(&allVars2,normSet,rangeName,code) ;
00186   
00187   return code+1 ; 
00188 } 
00189 
00190 
00191 
00192 //_____________________________________________________________________________
00193 Double_t RooProjectedPdf::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet*/, const char* rangeName) const 
00194 { 
00195   // Return analytical integral represent by appropriate element of projection cache
00196   
00197   CacheElem *cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
00198   
00199   if (cache) {
00200     Double_t ret= cache->_projection->getVal() ;
00201     return ret ;
00202   } else {
00203     
00204     RooArgSet* vars = getParameters(RooArgSet()) ;
00205     vars->add(intobs) ;
00206     RooArgSet* iset = _cacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
00207     RooArgSet* nset = _cacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
00208     
00209     Int_t code2(-1) ;
00210     const RooAbsReal* proj = getProjection(iset,nset,rangeName,code2) ;
00211     
00212     delete vars ;
00213     delete nset ;
00214     delete iset ;
00215     
00216     Double_t ret =  proj->getVal() ;
00217     return ret ;
00218   } 
00219   
00220 } 
00221 
00222 
00223 
00224 //_____________________________________________________________________________
00225 Int_t RooProjectedPdf::getGenerator(const RooArgSet& /*directVars*/, RooArgSet& /*generateVars*/, Bool_t /*staticInitOK*/) const 
00226  { 
00227    // No internal generator is implemented
00228    return 0 ; 
00229  } 
00230 
00231 
00232 
00233 //_____________________________________________________________________________
00234 void RooProjectedPdf::generateEvent(Int_t /*code*/) 
00235  { 
00236    // No internal generator is implemented
00237    return; 
00238  } 
00239 
00240 
00241 
00242 //_____________________________________________________________________________
00243 Bool_t RooProjectedPdf::redirectServersHook(const RooAbsCollection& newServerList, Bool_t /*mustReplaceAll*/, 
00244                                        Bool_t /*nameChange*/, Bool_t /*isRecursive*/) 
00245 {
00246   // Intercept a server redirection all and update list of dependents if necessary 
00247   // Specifically update the set proxy 'deps' which introduces the dependency
00248   // on server value dirty flags of ourselves
00249 
00250   // Redetermine explicit list of dependents if intPdf is being replaced
00251   RooAbsArg* newPdf = newServerList.find(intpdf.arg().GetName()) ;
00252   if (newPdf) {
00253 
00254     // Determine if set of dependens of new p.d.f is different from old p.d.f.
00255     RooArgSet olddeps(deps) ;
00256     RooArgSet* newdeps = newPdf->getParameters(intobs) ;
00257     RooArgSet* common = (RooArgSet*) newdeps->selectCommon(deps) ;    
00258     newdeps->remove(*common,kTRUE,kTRUE) ;
00259     olddeps.remove(*common,kTRUE,kTRUE) ;
00260 
00261     // If so, adjust composition of deps Listproxy
00262     if (newdeps->getSize()>0) {
00263       deps.add(*newdeps) ;
00264     }
00265     if (olddeps.getSize()>0) {
00266       deps.remove(olddeps,kTRUE,kTRUE) ;
00267     }
00268 
00269     delete common ;
00270     delete newdeps ;
00271   }
00272 
00273   return kFALSE ;
00274 }
00275 
00276 
00277 
00278 //_____________________________________________________________________________
00279 RooArgList RooProjectedPdf::CacheElem::containedArgs(Action)
00280 {
00281   // Return RooAbsArg elements contained in projection cache element.
00282   RooArgList ret(*_projection) ;  
00283   return ret ;
00284 }
00285 
00286 
00287 
00288 //_____________________________________________________________________________
00289 void RooProjectedPdf::printMetaArgs(ostream& os) const
00290 {
00291   // Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
00292   // integration operation
00293 
00294   os << "Int " << intpdf.arg().GetName() ;
00295  
00296   os << " d" ;
00297   os << intobs ;
00298   os << " " ;
00299   
00300 }
00301 
00302 
00303 
00304 
00305 //_____________________________________________________________________________
00306 void RooProjectedPdf::CacheElem::printCompactTreeHook(ostream& os, const char* indent, Int_t curElem, Int_t maxElem) 
00307 {
00308   // Print contents of cache when printing self as part of object tree
00309   
00310   if (curElem==0) {
00311     os << indent << "RooProjectedPdf begin projection cache" << endl ;
00312   }
00313 
00314   TString indent2(indent) ;
00315   indent2 += Form("[%d] ",curElem) ;
00316 
00317   _projection->printCompactTree(os,indent2) ;
00318 
00319   if(curElem==maxElem) {
00320     os << indent << "RooProjectedPdf end projection cache" << endl ;
00321   }
00322 }
00323 
00324 

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