Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooEffProd.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, NIKHEF
7  * GR, Gerhard Raven, NIKHEF/VU *
8  * *
9  * Redistribution and use in source and binary forms, *
10  * with or without modification, are permitted according to the terms *
11  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
12  *****************************************************************************/
13 
14 
15 /////////////////////////////////////////////////////////////////////////////////////
16 /// \class RooEffProd
17 /// The class RooEffProd implements the product of a PDF with an efficiency function.
18 /// The normalization integral of the product is calculated numerically, but the
19 /// event generation is handled by a specialized generator context that implements
20 /// the event generation in a more efficient for cases where the PDF has an internal
21 /// generator that is smarter than accept reject.
22 ///
23 
24 #include "RooFit.h"
25 #include "RooEffProd.h"
26 #include "RooEffGenContext.h"
27 #include "RooNameReg.h"
28 #include "RooRealVar.h"
29 
30 using namespace std;
31 
32 ClassImp(RooEffProd);
33  ;
34 
35 
36 
37 ////////////////////////////////////////////////////////////////////////////////
38 /// Constructor of a a production of p.d.f inPdf with efficiency
39 /// function inEff.
40 
41 RooEffProd::RooEffProd(const char *name, const char *title,
42  RooAbsPdf& inPdf, RooAbsReal& inEff) :
43  RooAbsPdf(name,title),
44  _cacheMgr(this,10),
45  _pdf("pdf","pre-efficiency pdf", this,inPdf),
46  _eff("eff","efficiency function",this,inEff),
47  _nset(0),
48  _fixedNset(0)
49 {
50 }
51 
52 
53 
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 /// Copy constructor
57 
58 RooEffProd::RooEffProd(const RooEffProd& other, const char* name) :
59  RooAbsPdf(other, name),
60  _cacheMgr(other._cacheMgr,this),
61  _pdf("pdf",this,other._pdf),
62  _eff("acc",this,other._eff),
63  _nset(0),
64  _fixedNset(0)
65 {
66 }
67 
68 
69 
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Destructor
73 
74 RooEffProd::~RooEffProd()
75 {
76 }
77 
78 
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Return p.d.f. value normalized over given set of observables
82 
83 Double_t RooEffProd::getValV(const RooArgSet* set) const
84 {
85  _nset = _fixedNset ? _fixedNset : set ;
86  return RooAbsPdf::getValV(set) ;
87 }
88 
89 
90 
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 /// Calculate and return 'raw' unnormalized value of p.d.f
94 
95 Double_t RooEffProd::evaluate() const
96 {
97  return eff()->getVal() * pdf()->getVal(_nset);
98 }
99 
100 
101 
102 ////////////////////////////////////////////////////////////////////////////////
103 /// Return specialized generator context for RooEffProds that implements generation
104 /// in a more efficient way than can be done for generic correlated products
105 
106 RooAbsGenContext* RooEffProd::genContext(const RooArgSet &vars, const RooDataSet *prototype,
107  const RooArgSet* auxProto, Bool_t verbose) const
108 {
109  assert(pdf()!=0);
110  assert(eff()!=0);
111  return new RooEffGenContext(*this,*pdf(),*eff(),vars,prototype,auxProto,verbose) ;
112 }
113 
114 
115 
116 
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// Return internal integration capabilities of the p.d.f. Given a set 'allVars' for which
120 /// integration is requested, returned the largest subset for which internal (analytical)
121 /// integration is implemented (in argument analVars). The return value is a unique integer
122 /// code that identifies the integration configuration (integrated observables and range name).
123 ///
124 /// This implementation in RooEffProd catches all integrals without normalization and reroutes them
125 /// through a custom integration routine that properly accounts for the use of normalized p.d.f.
126 /// in the evaluate() expression, which breaks the default RooAbsPdf normalization handling
127 
128 Int_t RooEffProd::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars,
129  const RooArgSet* normSet, const char* rangeName) const
130 {
131 
132  // No special handling required if a normalization set is given
133  if (normSet && normSet->getSize()>0) {
134  return 0 ;
135  }
136  // No special handling required if running with a fixed normalization set
137  if (_fixedNset) {
138  return 0 ;
139  }
140 
141  // Special handling of integrals w/o normalization set. We need to pass _a_ normalization set
142  // to pdf().getVal(nset) in evaluate() because for certain p.d.fs the shape depends on the
143  // chosen normalization set. This functions correctly automatically for plain getVal() calls,
144  // however when the (numeric) normalization is calculated, getVal() is called without any normalization
145  // set causing the normalization to be calculated for a possibly different shape. To fix that
146  // integrals over a RooEffProd without normalization setup are explicitly handled here. These integrals
147  // are calculated using a clone of the RooEffProd that has a fixed normalization set passed to the
148  // underlying p.d.f regardless of the normalization set passed to getVal(). Here the set of observables
149  // over which is integrated is passed.
150 
151  // Declare that we can analytically integrate all requested observables
152  analVars.add(allVars) ;
153 
154  // Check if cache was previously created
155  Int_t sterileIndex(-1) ;
156  CacheElem* cache = (CacheElem*) _cacheMgr.getObj(&allVars,&allVars,&sterileIndex,RooNameReg::ptr(rangeName)) ;
157  if (cache) {
158  return _cacheMgr.lastIndex()+1;
159  }
160 
161  // Construct cache with clone of p.d.f that has fixed normalization set that is passed to input pdf
162  cache = new CacheElem ;
163  cache->_intObs.addClone(allVars) ;
164  cache->_clone = (RooEffProd*) clone(Form("%s_clone",GetName())) ;
165  cache->_clone->_fixedNset = &cache->_intObs ;
166  cache->_int = cache->_clone->createIntegral(cache->_intObs,rangeName) ;
167 
168  // Store cache and return index as code
169  Int_t code = _cacheMgr.setObj(&allVars,&allVars,(RooAbsCacheElement*)cache,RooNameReg::ptr(rangeName)) ;
170 
171  return code+1 ;
172 }
173 
174 
175 
176 
177 
178 ////////////////////////////////////////////////////////////////////////////////
179 /// Return value of integral identified by code, which should be a return value of getAnalyticalIntegralWN,
180 /// Code zero is always handled and signifies no integration (return value is normalized p.d.f. value)
181 
182 Double_t RooEffProd::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
183 {
184  // Return analytical integral defined by given scenario code
185 
186  // No integration scenario
187  if (code==0) {
188  return getVal(normSet) ;
189  }
190 
191  // Partial integration scenarios
192  CacheElem* cache = (CacheElem*) _cacheMgr.getObjByIndex(code-1) ;
193 
194  return cache->_int->getVal() ;
195 }
196 
197 
198 
199 
200 ////////////////////////////////////////////////////////////////////////////////
201 /// Report all RooAbsArg derived objects contained in Cache Element (used in function optimization and
202 /// and server redirect management of the cache)
203 
204 RooArgList RooEffProd::CacheElem::containedArgs(Action)
205 {
206  RooArgList ret(_intObs) ;
207  ret.add(*_int) ;
208  ret.add(*_clone) ;
209  return ret ;
210 }