RooEffProd.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooEffProd.cxx 25184 2008-08-20 13:59:55Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, NIKHEF
00007  *   GR, Gerhard Raven, NIKHEF/VU                                            *
00008  *                                                                           *
00009  * Redistribution and use in source and binary forms,                        *
00010  * with or without modification, are permitted according to the terms        *
00011  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00012  *****************************************************************************/
00013 
00014 
00015 /////////////////////////////////////////////////////////////////////////////////////
00016 // BEGIN_HTML
00017 // The class RooEffProd implements the product of a PDF with an efficiency function.
00018 // The normalization integral of the product is calculated numerically, but the
00019 // event generation is handled by a specialized generator context that implements
00020 // the event generation in a more efficient for cases where the PDF has an internal
00021 // generator that is smarter than accept reject. 
00022 // END_HTML
00023 //
00024 
00025 #include "RooFit.h"
00026 #include "RooEffProd.h"
00027 #include "RooEffGenContext.h"
00028 #include "RooNameReg.h"
00029 #include "RooRealVar.h"
00030 
00031 ClassImp(RooEffProd)
00032   ;
00033 
00034 
00035 
00036 //_____________________________________________________________________________
00037 RooEffProd::RooEffProd(const char *name, const char *title, 
00038                              RooAbsPdf& inPdf, RooAbsReal& inEff) :
00039   RooAbsPdf(name,title),
00040   _cacheMgr(this,10),
00041   _pdf("pdf","pre-efficiency pdf", this,inPdf),
00042   _eff("eff","efficiency function",this,inEff),
00043   _nset(0),
00044   _fixedNset(0)
00045 {  
00046   // Constructor of a a production of p.d.f inPdf with efficiency
00047   // function inEff.
00048 }
00049 
00050 
00051 
00052 
00053 //_____________________________________________________________________________
00054 RooEffProd::RooEffProd(const RooEffProd& other, const char* name) : 
00055   RooAbsPdf(other, name),
00056   _cacheMgr(other._cacheMgr,this),
00057   _pdf("pdf",this,other._pdf),
00058   _eff("acc",this,other._eff),
00059   _nset(0),
00060   _fixedNset(0) 
00061 {
00062   // Copy constructor
00063 }
00064 
00065 
00066 
00067 
00068 //_____________________________________________________________________________
00069 RooEffProd::~RooEffProd() 
00070 {
00071   // Destructor
00072 }
00073 
00074 
00075 
00076 //_____________________________________________________________________________
00077 Double_t RooEffProd::getVal(const RooArgSet* set) const 
00078 {  
00079   // Return p.d.f. value normalized over given set of observables
00080 
00081   _nset = _fixedNset ? _fixedNset : set ;
00082   return RooAbsPdf::getVal(set) ;
00083 }
00084 
00085 
00086 
00087 
00088 //_____________________________________________________________________________
00089 Double_t RooEffProd::evaluate() const
00090 {
00091   // Calculate and return 'raw' unnormalized value of p.d.f
00092 
00093   return eff()->getVal() * pdf()->getVal(_nset);
00094 }
00095 
00096 
00097 
00098 //_____________________________________________________________________________
00099 RooAbsGenContext* RooEffProd::genContext(const RooArgSet &vars, const RooDataSet *prototype,
00100                                             const RooArgSet* auxProto, Bool_t verbose) const
00101 {
00102   // Return specialized generator context for RooEffProds that implements generation
00103   // in a more efficient way than can be done for generic correlated products
00104 
00105   assert(pdf()!=0);
00106   assert(eff()!=0);
00107   return new RooEffGenContext(*this,*pdf(),*eff(),vars,prototype,auxProto,verbose) ;
00108 }
00109 
00110 
00111 
00112 
00113 
00114 //_____________________________________________________________________________
00115 Int_t RooEffProd::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
00116                                           const RooArgSet* normSet, const char* rangeName) const 
00117 {
00118   // Return internal integration capabilities of the p.d.f. Given a set 'allVars' for which
00119   // integration is requested, returned the largest subset for which internal (analytical)
00120   // integration is implemented (in argument analVars). The return value is a unique integer
00121   // code that identifies the integration configuration (integrated observables and range name).
00122   //
00123   // This implementation in RooEffProd catches all integrals without normalization and reroutes them
00124   // through a custom integration routine that properly accounts for the use of normalized p.d.f.
00125   // in the evaluate() expression, which breaks the default RooAbsPdf normalization handling
00126 
00127   
00128   // No special handling required if a normalization set is given
00129   if (normSet && normSet->getSize()>0) {    
00130     return 0 ;
00131   }
00132   // No special handling required if running with a fixed normalization set
00133   if (_fixedNset) {
00134     return 0 ;
00135   }
00136 
00137   // Special handling of integrals w/o normalization set. We need to pass _a_ normalization set
00138   // to pdf().getVal(nset) in evaluate() because for certain p.d.fs the shape depends on the
00139   // chosen normalization set. This functions correctly automatically for plain getVal() calls,
00140   // however when the (numeric) normalization is calculated, getVal() is called without any normalization
00141   // set causing the normalization to be calculated for a possibly different shape. To fix that
00142   // integrals over a RooEffProd without normalization setup are explicitly handled here. These integrals
00143   // are calculated using a clone of the RooEffProd that has a fixed normalization set passed to the
00144   // underlying p.d.f regardless of the normalization set passed to getVal(). Here the set of observables
00145   // over which is integrated is passed.
00146 
00147   // Declare that we can analytically integrate all requested observables
00148   analVars.add(allVars) ;
00149 
00150   // Check if cache was previously created
00151   Int_t sterileIndex(-1) ;
00152   CacheElem* cache = (CacheElem*) _cacheMgr.getObj(&allVars,&allVars,&sterileIndex,RooNameReg::ptr(rangeName)) ;
00153   if (cache) {
00154     return _cacheMgr.lastIndex()+1;
00155   }
00156 
00157   // Construct cache with clone of p.d.f that has fixed normalization set that is passed to input pdf
00158   cache = new CacheElem ;
00159   cache->_intObs.addClone(allVars) ;
00160   cache->_clone = (RooEffProd*) clone(Form("%s_clone",GetName())) ;
00161   cache->_clone->_fixedNset = &cache->_intObs ;
00162   cache->_int = cache->_clone->createIntegral(cache->_intObs) ;
00163 
00164   // Store cache and return index as code
00165   Int_t code = _cacheMgr.setObj(&allVars,&allVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ; 
00166 
00167   return code+1 ;
00168 }
00169 
00170 
00171 
00172 
00173 
00174 //_____________________________________________________________________________
00175 Double_t RooEffProd::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const 
00176 {
00177   // Return value of integral identified by code, which should be a return value of getAnalyticalIntegralWN,
00178   // Code zero is always handled and signifies no integration (return value is normalized p.d.f. value)
00179 
00180   // Return analytical integral defined by given scenario code
00181 
00182   // No integration scenario
00183   if (code==0) {
00184     return getVal(normSet) ;
00185   }
00186 
00187   // Partial integration scenarios
00188   CacheElem* cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
00189 
00190   return cache->_int->getVal() ;
00191 }
00192 
00193 
00194 
00195 
00196 //_____________________________________________________________________________
00197 RooArgList RooEffProd::CacheElem::containedArgs(Action) 
00198 {
00199   // Report all RooAbsArg derived objects contained in Cache Element (used in function optimization and
00200   // and server redirect management of the cache)
00201 
00202   RooArgList ret(_intObs) ;
00203   ret.add(*_int) ;
00204   ret.add(*_clone) ;
00205   return ret ;
00206 }

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