Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooProjectedPdf.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2 * Project: RooFit *
3 * *
4 * Copyright (c) 2000-2005, Regents of the University of California *
5 * and Stanford University. All rights reserved. *
6 * *
7 * Redistribution and use in source and binary forms, *
8 * with or without modification, are permitted according to the terms *
9 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
10 *****************************************************************************/
11
12/**
13\file RooProjectedPdf.cxx
14\class RooProjectedPdf
15\ingroup Roofitcore
16
17A RooAbsPdf implementation that represent a projection
18of a given input p.d.f and the object returned by RooAbsPdf::createProjection.
19The actual projection integral for it value and normalization are
20calculated on the fly in getVal() once the normalization observables are known.
21Class RooProjectedPdf can cache projected p.d.f.s for multiple normalization
22observables simultaneously.
23The createProjection() method of RooProjectedPdf is overloaded and will
24return a new RooProjectedPdf that performs the projection of itself
25and the requested additional projections in one integration step
26The performance of <pre>f->createProjection(x)->createProjection(y)</pre>
27is therefore identical to that of <pre>f->createProjection(RooArgSet(x,y))</pre>
28**/
29
30#include "Riostream.h"
31
32#include "RooProjectedPdf.h"
33#include "RooMsgService.h"
34#include "RooAbsReal.h"
35#include "RooRealVar.h"
36#include "RooNameReg.h"
37
38
39
40////////////////////////////////////////////////////////////////////////////////
41/// Default constructor
42
44{
45}
46
47
48
49////////////////////////////////////////////////////////////////////////////////
50/// Construct projection of input pdf '_intpdf' over observables 'intObs'
51
52 RooProjectedPdf::RooProjectedPdf(const char *name, const char *title, RooAbsReal& _intpdf, const RooArgSet& intObs) :
53 RooAbsPdf(name,title),
54 intpdf("!IntegratedPdf","intpdf",this,_intpdf,false,false),
55 intobs("!IntegrationObservables","intobs",this,false,true),
56 deps("!Dependents","deps",this,true,false),
57 _cacheMgr(this,10)
58 {
59 // Since a projected PDF is an integral, we can use the same logic from
60 // RooRealIntegral via the projection integral to figure out what the
61 // servers are. Integration observables will be shape servers, the other
62 // servers are value servers.
63 int code;
64 auto proj = getProjection(&intObs, nullptr, nullptr, code);
65 for(RooAbsArg* server : proj->servers()) {
66 if(server->isShapeServer(*proj)) {
68 } else if(server->isValueServer(*proj)) {
69 deps.add(*server);
70 }
71 }
72 }
73
74
75
76////////////////////////////////////////////////////////////////////////////////
77/// Copy constructor
78
81 intpdf("!IntegratedPdf",this,other.intpdf),
82 intobs("!IntegrationObservable",this,other.intobs),
83 deps("!Dependents",this,other.deps),
84 _cacheMgr(other._cacheMgr,this)
85{
86 }
87
88
89
90////////////////////////////////////////////////////////////////////////////////
91/// Evaluate projected p.d.f
92
94{
95 // Calculate current unnormalized value of object
96 int code ;
97 return getProjection(&intobs, _normSet, normRange(), code)->getVal() ;
98}
99
100
101
102////////////////////////////////////////////////////////////////////////////////
103/// Retrieve object representing projection integral of input p.d.f
104/// over observables iset, while normalizing over observables
105/// nset. The code argument returned by reference is the unique code
106/// defining this particular projection configuration
107
108const RooAbsReal* RooProjectedPdf::getProjection(const RooArgSet* iset, const RooArgSet* nset, const char* rangeName, int& code) const
109{
110
111 // Check if this configuration was created before
112 Int_t sterileIdx(-1) ;
113 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObj(iset,nset,&sterileIdx,RooNameReg::ptr(rangeName)))) {
114 code = _cacheMgr.lastIndex() ;
115 return static_cast<const RooAbsReal*>(cache->_projection.get());
116 }
117
119 intpdf.arg().getObservables(nset, nset2);
120
121 if (iset) {
122 nset2.add(*iset) ;
123 }
124
125 auto cache = new CacheElem ;
126 cache->_projection = std::unique_ptr<RooAbsReal>{intpdf.arg().createIntegral(iset?*iset:RooArgSet(),&nset2,nullptr,rangeName)};
127
128 code = _cacheMgr.setObj(iset, nset, cache, RooNameReg::ptr(rangeName)) ;
129
130 coutI(Integration) << "RooProjectedPdf::getProjection(" << GetName() << ") creating new projection "
131 << cache->_projection->GetName() << " with code " << code << std::endl;
132
133 return cache->_projection.get();
134}
135
136
137
138////////////////////////////////////////////////////////////////////////////////
139/// Special version of RooAbsReal::createProjection that deals with
140/// projections of projections. Instead of integrating twice, a new
141/// RooProjectedPdf is returned that is configured to perform the
142/// complete integration in one step
143
145{
147 combiset.add(intobs) ;
148 return static_cast<RooAbsPdf&>( const_cast<RooAbsReal&>(intpdf.arg()) ).createProjection(combiset) ;
149}
150
151
152
153////////////////////////////////////////////////////////////////////////////////
154/// Force RooRealIntegral to relegate integration of all observables to internal logic
155
157{
158 return true ;
159}
160
161
162
163////////////////////////////////////////////////////////////////////////////////
164/// Mark all requested variables as internally integrated
165
167{
168 analVars.add(allVars) ;
169
170 // Create the appropriate integral
171 int code ;
172 RooArgSet allVars2(allVars) ;
173 allVars2.add(intobs) ;
175
176 return code+1 ;
177}
178
179
180
181////////////////////////////////////////////////////////////////////////////////
182/// Return analytical integral represent by appropriate element of projection cache
183
184double RooProjectedPdf::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet*/, const char* rangeName) const
185{
186 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObjByIndex(code-1))) {
187 return cache->_projection->getVal() ;
188 } else {
189
190 RooArgSet vars;
191 getParameters(nullptr, vars);
192 vars.add(intobs) ;
193 RooArgSet iset = _cacheMgr.selectFromSet1(vars, code-1) ;
194 RooArgSet nset = _cacheMgr.selectFromSet2(vars, code-1) ;
195
196 int code2 = -1 ;
197
198 return getProjection(&iset,&nset,rangeName,code2)->getVal() ;
199 }
200
201}
202
203
204
205////////////////////////////////////////////////////////////////////////////////
206/// Intercept a server redirection all and update list of dependents if necessary
207/// Specifically update the set proxy 'deps' which introduces the dependency
208/// on server value dirty flags of ourselves
209
211 bool nameChange, bool isRecursive)
212{
213 // Redetermine explicit list of dependents if intPdf is being replaced
214 if (RooAbsArg* newPdf = newServerList.find(intpdf.arg().GetName())) {
215
216 // Determine if set of dependents of new p.d.f is different from old p.d.f.
219 newPdf->getParameters(&intobs, newdeps);
221 newdeps.selectCommon(deps, common);
222 newdeps.remove(common,true,true) ;
223 olddeps.remove(common,true,true) ;
224
225 // If so, adjust composition of deps Listproxy
226 if (!newdeps.empty()) {
227 deps.add(newdeps) ;
228 }
229 if (!olddeps.empty()) {
230 deps.remove(olddeps,true,true) ;
231 }
232 }
233
235}
236
237
238
239////////////////////////////////////////////////////////////////////////////////
240/// Return RooAbsArg elements contained in projection cache element.
241
246
247
248
249////////////////////////////////////////////////////////////////////////////////
250/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
251/// integration operation
252
253void RooProjectedPdf::printMetaArgs(std::ostream& os) const
254{
255 os << "Int " << intpdf.arg().GetName() << " d" << intobs << " ";
256}
257
258
259
260
261////////////////////////////////////////////////////////////////////////////////
262/// Print contents of cache when printing self as part of object tree
263
265{
266 if (curElem==0) {
267 os << indent << "RooProjectedPdf begin projection cache" << std::endl ;
268 }
269
271 indent2 += Form("[%d] ",curElem) ;
272
273 _projection->printCompactTree(os,indent2) ;
274
275 if(curElem==maxElem) {
276 os << indent << "RooProjectedPdf end projection cache" << std::endl ;
277 }
278}
279
280
281std::unique_ptr<RooAbsArg>
283{
286 nset2.add(intobs);
287
288 auto newArg = std::unique_ptr<RooAbsReal>{intpdf->createIntegral(intobs, &nset2)};
290 return newArg;
291}
#define coutI(a)
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
Abstract container object that can hold multiple RooAbsArg objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:320
const char * normRange() const
Definition RooAbsPdf.h:250
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
RooArgSet selectFromSet1(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 1 with a given index and an i...
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet containing the objects that are both in the cached set 2 with a given index and an i...
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false) override
Remove object 'var' from set and deregister 'var' as server to owner.
void markAsCompiled(RooAbsArg &arg) const
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
RooArgList containedArgs(Action) override
Return RooAbsArg elements contained in projection cache element.
void printCompactTreeHook(std::ostream &, const char *, Int_t, Int_t) override
Print contents of cache when printing self as part of object tree.
std::unique_ptr< RooAbsReal > _projection
A RooAbsPdf implementation that represent a projection of a given input p.d.f and the object returned...
RooRealProxy intpdf
p.d.f that is integrated
RooObjCacheManager _cacheMgr
! The cache manager
double evaluate() const override
Evaluate projected p.d.f.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integral represent by appropriate element of projection cache.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
Force RooRealIntegral to relegate integration of all observables to internal logic.
bool redirectServersHook(const RooAbsCollection &newServerList, bool, bool, bool) override
Intercept a server redirection all and update list of dependents if necessary Specifically update the...
RooAbsPdf * createProjection(const RooArgSet &iset) override
Special version of RooAbsReal::createProjection that deals with projections of projections.
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the...
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
RooProjectedPdf()
Default constructor.
const RooAbsReal * getProjection(const RooArgSet *iset, const RooArgSet *nset, const char *rangeName, int &code) const
Retrieve object representing projection integral of input p.d.f over observables iset,...
RooSetProxy intobs
observables that p.d.f is integrated over
RooSetProxy deps
dependents of this p.d.f
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Mark all requested variables as internally integrated.
const T & arg() const
Return reference to object held in proxy.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Basic string class.
Definition TString.h:139