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#include "RooWrapperPdf.h"
38
39 ////////////////////////////////////////////////////////////////////////////////
40 /// Default constructor
41
43
44 ////////////////////////////////////////////////////////////////////////////////
45 /// Construct projection of input pdf '_intpdf' over observables 'intObs'
46
47 RooProjectedPdf::RooProjectedPdf(const char *name, const char *title, RooAbsReal& _intpdf, const RooArgSet& intObs) :
48 RooAbsPdf(name,title),
49 intpdf("!IntegratedPdf","intpdf",this,_intpdf,false,false),
50 intobs("!IntegrationObservables","intobs",this,false,true),
51 deps("!Dependents","deps",this,true,false),
52 _cacheMgr(this,10)
53 {
54 // Since a projected PDF is an integral, we can use the same logic from
55 // RooRealIntegral via the projection integral to figure out what the
56 // servers are. Integration observables will be shape servers, the other
57 // servers are value servers.
58 int code;
59 auto proj = getProjection(&intObs, nullptr, nullptr, code);
60 for(RooAbsArg* server : proj->servers()) {
61 if(server->isShapeServer(*proj)) {
63 } else if(server->isValueServer(*proj)) {
64 deps.add(*server);
65 }
66 }
67 }
68
69
70
71////////////////////////////////////////////////////////////////////////////////
72/// Copy constructor
73
76 intpdf("!IntegratedPdf",this,other.intpdf),
77 intobs("!IntegrationObservables",this,other.intobs),
78 deps("!Dependents",this,other.deps),
79 _cacheMgr(other._cacheMgr,this)
80{
81 }
82
83
84
85////////////////////////////////////////////////////////////////////////////////
86/// Evaluate projected p.d.f
87
89{
90 // Calculate current unnormalized value of object
91 int code ;
92 return getProjection(&intobs, _normSet, normRange(), code)->getVal() ;
93}
94
95
96
97////////////////////////////////////////////////////////////////////////////////
98/// Retrieve object representing projection integral of input p.d.f
99/// over observables iset, while normalizing over observables
100/// nset. The code argument returned by reference is the unique code
101/// defining this particular projection configuration
102
103const RooAbsReal* RooProjectedPdf::getProjection(const RooArgSet* iset, const RooArgSet* nset, const char* rangeName, int& code) const
104{
105
106 // Check if this configuration was created before
107 Int_t sterileIdx(-1) ;
108 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObj(iset,nset,&sterileIdx,RooNameReg::ptr(rangeName)))) {
109 code = _cacheMgr.lastIndex() ;
110 return static_cast<const RooAbsReal*>(cache->_projection.get());
111 }
112
114 intpdf.arg().getObservables(nset, nset2);
115
116 if (iset) {
117 nset2.add(*iset) ;
118 }
119
120 auto cache = new CacheElem ;
121 cache->_projection = std::unique_ptr<RooAbsReal>{intpdf.arg().createIntegral(iset?*iset:RooArgSet(),&nset2,nullptr,rangeName)};
122
123 code = _cacheMgr.setObj(iset, nset, cache, RooNameReg::ptr(rangeName)) ;
124
125 coutI(Integration) << "RooProjectedPdf::getProjection(" << GetName() << ") creating new projection "
126 << cache->_projection->GetName() << " with code " << code << std::endl;
127
128 return cache->_projection.get();
129}
130
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Special version of RooAbsReal::createProjection that deals with
135/// projections of projections. Instead of integrating twice, a new
136/// RooProjectedPdf is returned that is configured to perform the
137/// complete integration in one step
138
140{
142 combiset.add(intobs) ;
143 return static_cast<RooAbsPdf&>( const_cast<RooAbsReal&>(intpdf.arg()) ).createProjection(combiset) ;
144}
145
146
147
148////////////////////////////////////////////////////////////////////////////////
149/// Force RooRealIntegral to relegate integration of all observables to internal logic
150
152{
153 return true ;
154}
155
156
157
158////////////////////////////////////////////////////////////////////////////////
159/// Mark all requested variables as internally integrated
160
162{
163 analVars.add(allVars) ;
164
165 // Create the appropriate integral
166 int code ;
167 RooArgSet allVars2(allVars) ;
168 allVars2.add(intobs) ;
170
171 return code+1 ;
172}
173
174
175
176////////////////////////////////////////////////////////////////////////////////
177/// Return analytical integral represent by appropriate element of projection cache
178
179double RooProjectedPdf::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet*/, const char* rangeName) const
180{
181 if (auto cache = static_cast<CacheElem*>(_cacheMgr.getObjByIndex(code-1))) {
182 return cache->_projection->getVal() ;
183 } else {
184
185 RooArgSet vars;
186 getParameters(nullptr, vars);
187 vars.add(intobs) ;
188 RooArgSet iset = _cacheMgr.selectFromSet1(vars, code-1) ;
189 RooArgSet nset = _cacheMgr.selectFromSet2(vars, code-1) ;
190
191 int code2 = -1 ;
192
193 return getProjection(&iset,&nset,rangeName,code2)->getVal() ;
194 }
195
196}
197
198
199
200////////////////////////////////////////////////////////////////////////////////
201/// Intercept a server redirection all and update list of dependents if necessary
202/// Specifically update the set proxy 'deps' which introduces the dependency
203/// on server value dirty flags of ourselves
204
206 bool nameChange, bool isRecursive)
207{
208 // Redetermine explicit list of dependents if intPdf is being replaced
209 if (RooAbsArg* newPdf = newServerList.find(intpdf.arg().GetName())) {
210
211 // Determine if set of dependents of new p.d.f is different from old p.d.f.
214 newPdf->getParameters(&intobs, newdeps);
216 newdeps.selectCommon(deps, common);
217 newdeps.remove(common,true,true) ;
218 olddeps.remove(common,true,true) ;
219
220 // If so, adjust composition of deps Listproxy
221 if (!newdeps.empty()) {
222 deps.add(newdeps) ;
223 }
224 if (!olddeps.empty()) {
225 deps.remove(olddeps,true,true) ;
226 }
227 }
228
230}
231
232
233
234////////////////////////////////////////////////////////////////////////////////
235/// Return RooAbsArg elements contained in projection cache element.
236
241
242
243
244////////////////////////////////////////////////////////////////////////////////
245/// Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
246/// integration operation
247
248void RooProjectedPdf::printMetaArgs(std::ostream& os) const
249{
250 os << "Int " << intpdf.arg().GetName() << " d" << intobs << " ";
251}
252
253
254
255
256////////////////////////////////////////////////////////////////////////////////
257/// Print contents of cache when printing self as part of object tree
258
260{
261 if (curElem==0) {
262 os << indent << "RooProjectedPdf begin projection cache" << std::endl ;
263 }
264
266 indent2 += Form("[%d] ",curElem) ;
267
268 _projection->printCompactTree(os,indent2) ;
269
270 if(curElem==maxElem) {
271 os << indent << "RooProjectedPdf end projection cache" << std::endl ;
272 }
273}
274
275
276std::unique_ptr<RooAbsArg>
278{
281 nset2.add(intobs);
282
283 auto newArg = std::unique_ptr<RooAbsReal>{intpdf->createIntegral(intobs, &nset2)};
284
285 std::string namePdf = std::string{newArg->GetName()} + "_wrapped_pdf";
286
287 auto newArgPdf = std::make_unique<RooWrapperPdf>(namePdf.c_str(), namePdf.c_str(), *newArg);
288
291
292 newArgPdf->addOwnedComponents(std::move(newArg));
293
294 return newArgPdf;
295}
#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:2496
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
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:32
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:316
const char * normRange() const
Definition RooAbsPdf.h:246
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:63
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:138