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
17Class RooProjectedPdf is a 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
39
40
41////////////////////////////////////////////////////////////////////////////////
42/// Default constructor
43
45{
46}
47
48
49
50////////////////////////////////////////////////////////////////////////////////
51/// Construct projection of input pdf '_intpdf' over observables 'intObs'
52
53 RooProjectedPdf::RooProjectedPdf(const char *name, const char *title, RooAbsReal& _intpdf, const RooArgSet& intObs) :
54 RooAbsPdf(name,title),
55 intpdf("!IntegratedPdf","intpdf",this,_intpdf,false,false),
56 intobs("!IntegrationObservables","intobs",this,false,true),
57 deps("!Dependents","deps",this,true,false),
58 _cacheMgr(this,10)
59 {
60 // Since a projected PDF is an integral, we can use the same logic from
61 // RooRealIntegral via the projection integral to figure out what the
62 // servers are. Integration observables will be shape servers, the other
63 // servers are value servers.
64 int code;
65 auto proj = getProjection(&intObs, nullptr, 0, code);
66 for(RooAbsArg* server : proj->servers()) {
67 if(server->isShapeServer(*proj)) {
68 intobs.add(*server);
69 } else if(server->isValueServer(*proj)) {
70 deps.add(*server);
71 }
72 }
73 }
74
75
76
77////////////////////////////////////////////////////////////////////////////////
78/// Copy constructor
79
81 RooAbsPdf(other,name),
82 intpdf("!IntegratedPdf",this,other.intpdf),
83 intobs("!IntegrationObservable",this,other.intobs),
84 deps("!Dependents",this,other.deps),
85 _cacheMgr(other._cacheMgr,this)
86{
87 }
88
89
90
91////////////////////////////////////////////////////////////////////////////////
92/// Evaluate projected p.d.f
93
95{
96 // Calculate current unnormalized value of object
97 int code ;
98 return getProjection(&intobs, _normSet, normRange(), code)->getVal() ;
99}
100
101
102
103////////////////////////////////////////////////////////////////////////////////
104/// Retrieve object representing projection integral of input p.d.f
105/// over observables iset, while normalizing over observables
106/// nset. The code argument returned by reference is the unique code
107/// defining this particular projection configuration
108
109const RooAbsReal* RooProjectedPdf::getProjection(const RooArgSet* iset, const RooArgSet* nset, const char* rangeName, int& code) const
110{
111
112 // Check if this configuration was created before
113 Int_t sterileIdx(-1) ;
114 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObj(iset,nset,&sterileIdx,RooNameReg::ptr(rangeName)))) {
115 code = _cacheMgr.lastIndex() ;
116 return static_cast<const RooAbsReal*>(cache->_projection.get());
117 }
118
119 RooArgSet nset2;
120 intpdf.arg().getObservables(nset, nset2);
121
122 if (iset) {
123 nset2.add(*iset) ;
124 }
125
126 auto cache = new CacheElem ;
127 cache->_projection = std::unique_ptr<RooAbsReal>{intpdf.arg().createIntegral(iset?*iset:RooArgSet(),&nset2,0,rangeName)};
128
129 code = _cacheMgr.setObj(iset, nset, cache, RooNameReg::ptr(rangeName)) ;
130
131 coutI(Integration) << "RooProjectedPdf::getProjection(" << GetName() << ") creating new projection "
132 << cache->_projection->GetName() << " with code " << code << std::endl;
133
134 return cache->_projection.get();
135}
136
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// Special version of RooAbsReal::createProjection that deals with
141/// projections of projections. Instead of integrating twice, a new
142/// RooProjectedPdf is returned that is configured to perform the
143/// complete integration in one step
144
146{
147 RooArgSet combiset(iset) ;
148 combiset.add(intobs) ;
149 return static_cast<RooAbsPdf&>( const_cast<RooAbsReal&>(intpdf.arg()) ).createProjection(combiset) ;
150}
151
152
153
154////////////////////////////////////////////////////////////////////////////////
155/// Force RooRealIntegral to relegate integration of all observables to internal logic
156
158{
159 return true ;
160}
161
162
163
164////////////////////////////////////////////////////////////////////////////////
165/// Mark all requested variables as internally integrated
166
167Int_t RooProjectedPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName) const
168{
169 analVars.add(allVars) ;
170
171 // Create the appropriate integral
172 int code ;
173 RooArgSet allVars2(allVars) ;
174 allVars2.add(intobs) ;
175 getProjection(&allVars2,normSet,rangeName,code) ;
176
177 return code+1 ;
178}
179
180
181
182////////////////////////////////////////////////////////////////////////////////
183/// Return analytical integral represent by appropriate element of projection cache
184
185double RooProjectedPdf::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet*/, const char* rangeName) const
186{
187 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObjByIndex(code-1))) {
188 return cache->_projection->getVal() ;
189 } else {
190
191 RooArgSet vars;
192 getParameters(nullptr, vars);
193 vars.add(intobs) ;
194 RooArgSet iset = _cacheMgr.selectFromSet1(vars, code-1) ;
195 RooArgSet nset = _cacheMgr.selectFromSet2(vars, code-1) ;
196
197 int code2 = -1 ;
198
199 return getProjection(&iset,&nset,rangeName,code2)->getVal() ;
200 }
201
202}
203
204
205
206////////////////////////////////////////////////////////////////////////////////
207/// Intercept a server redirection all and update list of dependents if necessary
208/// Specifically update the set proxy 'deps' which introduces the dependency
209/// on server value dirty flags of ourselves
210
211bool RooProjectedPdf::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll,
212 bool nameChange, bool isRecursive)
213{
214 // Redetermine explicit list of dependents if intPdf is being replaced
215 if (RooAbsArg* newPdf = newServerList.find(intpdf.arg().GetName())) {
216
217 // Determine if set of dependens of new p.d.f is different from old p.d.f.
218 RooArgSet olddeps(deps) ;
219 RooArgSet newdeps;
220 newPdf->getParameters(&intobs, newdeps);
221 RooArgSet common;
222 newdeps.selectCommon(deps, common);
223 newdeps.remove(common,true,true) ;
224 olddeps.remove(common,true,true) ;
225
226 // If so, adjust composition of deps Listproxy
227 if (!newdeps.empty()) {
228 deps.add(newdeps) ;
229 }
230 if (!olddeps.empty()) {
231 deps.remove(olddeps,true,true) ;
232 }
233 }
234
235 return RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
236}
237
238
239
240////////////////////////////////////////////////////////////////////////////////
241/// Return RooAbsArg elements contained in projection cache element.
242
244{
245 return *_projection;
246}
247
248
249
250////////////////////////////////////////////////////////////////////////////////
251/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
252/// integration operation
253
254void RooProjectedPdf::printMetaArgs(std::ostream& os) const
255{
256 os << "Int " << intpdf.arg().GetName() << " d" << intobs << " ";
257}
258
259
260
261
262////////////////////////////////////////////////////////////////////////////////
263/// Print contents of cache when printing self as part of object tree
264
265void RooProjectedPdf::CacheElem::printCompactTreeHook(std::ostream& os, const char* indent, Int_t curElem, Int_t maxElem)
266{
267 if (curElem==0) {
268 os << indent << "RooProjectedPdf begin projection cache" << std::endl ;
269 }
270
271 TString indent2(indent) ;
272 indent2 += Form("[%d] ",curElem) ;
273
274 _projection->printCompactTree(os,indent2) ;
275
276 if(curElem==maxElem) {
277 os << indent << "RooProjectedPdf end projection cache" << std::endl ;
278 }
279}
280
281
282std::unique_ptr<RooAbsArg>
284{
285 RooArgSet nset2;
286 intpdf->getObservables(&normSet, nset2);
287 nset2.add(intobs);
288
289 auto newArg = std::unique_ptr<RooAbsReal>{intpdf->createIntegral(intobs, &nset2)};
290 ctx.markAsCompiled(*newArg);
291 return newArg;
292}
#define coutI(a)
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
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.
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
bool selectCommon(const RooAbsCollection &refColl, RooAbsCollection &outColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:377
const char * normRange() const
Definition RooAbsPdf.h:310
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.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
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:55
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 contatining the objects that are both in the cached set 1.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
RooArgSet selectFromSet2(RooArgSet const &argSet, int index) const
Create RooArgSet contatining the objects that are both in the cached set 2.
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
Class RooProjectedPdf is a RooAbsPdf implementation that represent a projection of a given input p....
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:47
Basic string class.
Definition TString.h:139