Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooNumConvolution.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, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooNumConvolution.cxx
19\class RooNumConvolution
20\ingroup Roofitcore
21
22Numeric 1-dimensional convolution operator PDF. This class can convolve any PDF
23with any other PDF
24This class should not be used blindly as numeric convolution is computing
25intensive and prone to stability fitting problems. If an analytic convolution
26can be calculated, you should use that or implement it if not available.
27RooNumConvolution implements reasonable defaults that should convolve most
28functions reasonably well, but results strongly depend on the shape of your
29input PDFS so always check your result.
30
31The default integration engine for the numeric convolution is the
32adaptive Gauss-Kronrod method, which empirically seems the most robust
33for this task. You can override the convolution integration settings via
34the RooNumIntConfig object reference returned by the convIntConfig() member
35function
36By default the numeric convolution is integrated from -infinity to
37+infinity through a <pre>x -> 1/x</pre> coordinate transformation of the
38tails. For convolution with a very small bandwidth it may be
39advantageous (for both CPU consumption and stability) if the
40integration domain is limited to a finite range. The function
41setConvolutionWindow(mean,width,scale) allows to set a sliding
42window around the x value to be calculated taking a RooAbsReal
43expression for an offset and a width to be taken around the x
44value. These input expression can be RooFormulaVars or other
45function objects although the 3d 'scale' argument 'scale'
46multiplies the width RooAbsReal expression given in the 2nd
47argument, allowing for an appropriate window definition for most
48cases without need for a RooFormulaVar object: e.g. a Gaussian
49resolution PDF do setConvolutionWindow(gaussMean,gaussSigma,5)
50Note that for a 'wide' Gaussian the -inf to +inf integration
51may converge more quickly than that over a finite range!
52The default numeric precision is 1e-7, i.e. the global default for
53numeric integration but you should experiment with this value to
54see if it is sufficient for example by studying the number of function
55calls that MINUIT needs to fit your function as function of the
56convolution precision.
57**/
58
59#include "Riostream.h"
60#include "TH2F.h"
61#include "RooNumConvolution.h"
62#include "RooArgList.h"
63#include "RooRealVar.h"
64#include "RooCustomizer.h"
66#include "RooNumIntFactory.h"
67#include "RooGenContext.h"
68#include "RooConvGenContext.h"
69#include "RooMsgService.h"
70
71
72using std::endl, std::ostream;
73
75
76
77
78////////////////////////////////////////////////////////////////////////////////
79
81 _init(false),
82 _integrand(nullptr),
83 _cloneVar(nullptr),
84 _clonePdf(nullptr),
85 _cloneModel(nullptr),
86 _useWindow(false),
87 _windowScale(1),
88 _verboseThresh(2000),
89 _doProf(false),
90 _callHist(nullptr)
91{
92}
93
94
95
96////////////////////////////////////////////////////////////////////////////////
97/// Constructor of convolution operator PDF
98///
99/// convVar : convolution variable (on which both pdf and resmodel should depend)
100/// pdf : input 'physics' pdf
101/// resmodel : input 'resolution' pdf
102///
103/// output is pdf(x) (X) resmodel(x) = Int [ pdf(x') resmodel (x-x') ] dx'
104///
105
106RooNumConvolution::RooNumConvolution(const char *name, const char *title, RooRealVar& convVar, RooAbsReal& inPdf, RooAbsReal& resmodel, const RooNumConvolution* proto) :
107 RooAbsReal(name,title),
108 _init(false),
109 _convIntConfig(RooNumIntConfig::defaultConfig()),
110 _integrand(nullptr),
111 _origVar("origVar","Original Convolution variable",this,convVar),
112 _origPdf("origPdf","Original Input PDF",this,inPdf),
113 _origModel("origModel","Original Resolution model",this,resmodel),
114 _ownedClonedPdfSet("ownedClonePdfSet"),
115 _ownedClonedModelSet("ownedCloneModelSet"),
116 _cloneVar(nullptr),
117 _clonePdf(nullptr),
118 _cloneModel(nullptr),
119 _useWindow(false),
120 _windowScale(1),
121 _windowParam("windowParam","Convolution window parameter",this,false),
122 _verboseThresh(2000),
123 _doProf(false),
124 _callHist(nullptr)
125{
126 // Use Adaptive Gauss-Kronrod integration by default for the convolution integral
127 _convIntConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
128 _convIntConfig.method1DOpen().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
129
130 if (proto) {
131 convIntConfig() = proto->convIntConfig() ;
132 if (proto->_useWindow) {
133 setConvolutionWindow(static_cast<RooAbsReal&>(*proto->_windowParam.at(0)),static_cast<RooAbsReal&>(*proto->_windowParam.at(1)),proto->_windowScale) ;
134 }
135 }
136}
137
138
139
140////////////////////////////////////////////////////////////////////////////////
141/// Copy constructor
142
144 RooAbsReal(other,name),
145 _init(false),
146 _convIntConfig(other._convIntConfig),
147 _integrand(nullptr),
148 _origVar("origVar",this,other._origVar),
149 _origPdf("origPdf",this,other._origPdf),
150 _origModel("origModel",this,other._origModel),
151 _ownedClonedPdfSet("ownedClonePdfSet"),
152 _ownedClonedModelSet("ownedCloneModelSet"),
153 _cloneVar(nullptr),
154 _clonePdf(nullptr),
155 _cloneModel(nullptr),
156 _useWindow(other._useWindow),
157 _windowScale(other._windowScale),
158 _windowParam("windowParam",this,other._windowParam),
159 _verboseThresh(other._verboseThresh),
160 _doProf(other._doProf),
161 _callHist(other._callHist)
162{
163}
164
165
166
167////////////////////////////////////////////////////////////////////////////////
168/// One-time initialization of object
169
171{
172 // Initialization function -- create clone of convVar (x') and deep-copy clones of pdf and
173 // model that are connected to x' rather than x (convVar)
174
175 // Start out clean
178
179 if (_cloneVar) delete _cloneVar ;
180
181 // Customize a copy of origPdf that is connected to x' rather than x
182 // store all cloned components in _clonePdfSet as well as x' itself
183 _cloneVar = new RooRealVar(Form("%s_prime",_origVar.arg().GetName()),"Convolution Variable",0) ;
184
185 RooCustomizer mgr1(pdf(),"NumConv_PdfClone") ;
187 mgr1.replaceArg(var(),*_cloneVar) ;
188 _clonePdf = static_cast<RooAbsReal*>(mgr1.build()) ;
189
190 RooCustomizer mgr2(model(),"NumConv_ModelClone") ;
192 mgr2.replaceArg(var(),*_cloneVar) ;
193 _cloneModel = static_cast<RooAbsReal*>(mgr2.build()) ;
194
195 // Change name back to original name
197
198 // Create Convolution integrand
200
201 // Instantiate integrator for convolution integrand
203 _integrator->setUseIntegrandLimits(false) ;
204
205 _init = true ;
206}
207
208
209
210
211////////////////////////////////////////////////////////////////////////////////
212/// Destructor
213
215{
216}
217
218
219
220////////////////////////////////////////////////////////////////////////////////
221/// Calculate convolution integral
222
224{
225 // Check if deferred initialization has occurred
226 if (!_init) initialize() ;
227
228 // Retrieve current value of convolution variable
229 double x = _origVar ;
230
231 // Propagate current normalization set to integrand
233
234 // Adjust convolution integration window
235 if (_useWindow) {
236 double center = (static_cast<RooAbsReal*>(_windowParam.at(0)))->getVal() ;
237 double width = _windowScale * (static_cast<RooAbsReal*>(_windowParam.at(1)))->getVal() ;
238 _integrator->setLimits(x-center-width,x-center+width) ;
239 } else {
241 }
242
243 // Calculate convolution for present x
245 double ret = _integrator->integral(&x) ;
246 if (_doProf) {
249 coutW(Integration) << "RooNumConvolution::evaluate(" << GetName() << ") WARNING convolution integral at x=" << x
250 << " required " << _integrand->numCall() << " function evaluations" << endl ;
251 }
252 }
253
254 return ret ;
255}
256
257
258
259////////////////////////////////////////////////////////////////////////////////
260/// Intercept server redirects. Throw away cache, as figuring out redirections on the cache is an unsolvable problem.
261
262bool RooNumConvolution::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll,
263 bool nameChange, bool isRecursive)
264{
265 _init = false ;
266 return RooAbsReal::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
267}
268
269
270
271////////////////////////////////////////////////////////////////////////////////
272/// Removes previously defined convolution window, reverting to convolution from -inf to +inf
273
275{
276 _useWindow = false ;
278}
279
280
281
282////////////////////////////////////////////////////////////////////////////////
283/// Restrict convolution integral to finite range [ x - C - S*W, x - C + S*W ]
284/// where x is current value of convolution variablem, C = centerParam, W=widthParam and S = widthScaleFactor
285/// Inputs centerParam and withParam can be function expressions (RooAbsReal, RooFormulaVar) etc.
286
287void RooNumConvolution::setConvolutionWindow(RooAbsReal& centerParam, RooAbsReal& widthParam, double widthScaleFactor)
288{
289 _useWindow = true ;
291 _windowParam.add(centerParam) ;
292 _windowParam.add(widthParam) ;
293 _windowScale = widthScaleFactor ;
294}
295
296
297
298////////////////////////////////////////////////////////////////////////////////
299/// Activate warning messages if number of function calls needed for evaluation of convolution integral
300/// exceeds given threshold
301
303{
304 if (threshold<0) {
305 coutE(InputArguments) << "RooNumConvolution::setCallWarning(" << GetName() << ") ERROR: threshold must be positive, value unchanged" << endl ;
306 return ;
307 }
308 _verboseThresh = threshold ;
309}
310
311
312////////////////////////////////////////////////////////////////////////////////
313/// Activate call profile if flag is set to true. A 2-D histogram is kept that stores the required number
314/// of function calls versus the value of x, the convolution variable
315///
316/// All clones of RooNumConvolution objects will keep logging to the histogram of the original class
317/// so that performance of temporary object clones, such as used in e.g. fitting, plotting and generating
318/// are all logged in a single place.
319///
320/// Function caller should take ownership of profiling histogram as it is not deleted at the RooNumConvolution dtor
321///
322/// Calling this function with flag set to false will deactivate call profiling and delete the profiling histogram
323
324void RooNumConvolution::setCallProfiling(bool flag, Int_t nbinX, Int_t nbinCall, Int_t nCallHigh)
325{
326 if (flag) {
327 if (_doProf) {
328 delete _callHist ;
329 }
330 _callHist = new TH2F(Form("callHist_%s",GetName()),Form("Call Profiling of RooNumConvolution %s",GetTitle()),
331 nbinX,_origVar.min(),_origVar.max(),
332 nbinCall,0,nCallHigh) ;
333 _doProf=true ;
334
335 } else if (_doProf) {
336
337 delete _callHist ;
338 _callHist = nullptr ;
339 _doProf = false ;
340 }
341
342}
343
344
345
346////////////////////////////////////////////////////////////////////////////////
347/// Hook function to intercept printCompactTree() calls so that it can print out
348/// the content of its private cache in the print sequence
349
351{
352 os << indent << "RooNumConvolution begin cache" << endl ;
353
354 if (_init) {
355 _cloneVar->printCompactTree(os,Form("%s[Var]",indent)) ;
356 _clonePdf->printCompactTree(os,Form("%s[Pdf]",indent)) ;
357 _cloneModel->printCompactTree(os,Form("%s[Mod]",indent)) ;
358 }
359
360 os << indent << "RooNumConvolution end cache" << endl ;
361}
362
363
bool _init
! Is object initialized
#define coutW(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
Option_t Option_t width
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
const char * proto
Definition civetweb.c:17536
void SetName(const char *name) override
Set the name of the TNamed.
void printCompactTree(const char *indent="", const char *fileName=nullptr, const char *namePat=nullptr, RooAbsArg *client=nullptr)
Print tree structure of expression tree on stdout, or to file if filename is specified.
Abstract container object that can hold multiple RooAbsArg objects.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
void resetNumCall() const
Reset function call counter.
Definition RooAbsFunc.h:52
Int_t numCall() const
Return number of function calls since last reset.
Definition RooAbsFunc.h:47
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Function that is called at the end of redirectServers().
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
bool setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
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...
Implementation of RooAbsFunc that represent the integrand of a generic (numeric) convolution A (x) B ...
void setNormalizationSet(const RooArgSet *nset)
RooCustomizer is a factory class to produce clones of a prototype composite PDF object with the same ...
void setCloneBranchSet(RooArgSet &cloneBranchSet)
Releases ownership of list of cloned branch nodes.
void replaceArg(const RooAbsArg &orig, const RooAbsArg &subst)
Replace any occurrence of arg 'orig' with arg 'subst'.
RooAbsArg * build(const char *masterCatState, bool verbose=false)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
Numeric 1-dimensional convolution operator PDF.
void clearConvolutionWindow()
Removes previously defined convolution window, reverting to convolution from -inf to +inf.
RooArgSet _ownedClonedModelSet
Owning set of cloned model components.
std::unique_ptr< RooAbsIntegrator > _integrator
! Numeric integrator of convolution integrand
RooAbsReal * _cloneModel
Pointer to cloned model.
void setCallProfiling(bool flag, Int_t nbinX=40, Int_t nbinCall=40, Int_t nCallHigh=1000)
Activate call profile if flag is set to true.
~RooNumConvolution() override
Destructor.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override
Intercept server redirects. Throw away cache, as figuring out redirections on the cache is an unsolva...
void initialize() const
One-time initialization of object.
RooRealProxy _origVar
Original convolution variable.
Int_t _verboseThresh
Call count threshold for verbose printing.
RooConvIntegrandBinding * _integrand
! Binding of Convolution Integrand function
void printCompactTreeHook(std::ostream &os, const char *indent="") override
Hook function to intercept printCompactTree() calls so that it can print out the content of its priva...
double evaluate() const override
Calculate convolution integral.
RooListProxy _windowParam
Holder for optional convolution integration window scaling parameter.
double _windowScale
Scale factor for window parameter.
RooArgSet _ownedClonedPdfSet
Owning set of cloned PDF components.
RooRealVar & var() const
bool _doProf
Switch to activate profiling option.
RooAbsReal * _cloneVar
Pointer to cloned convolution variable.
RooAbsReal * _clonePdf
Pointer to cloned PDF.
RooNumIntConfig & convIntConfig()
void setCallWarning(Int_t threshold=2000)
Activate warning messages if number of function calls needed for evaluation of convolution integral e...
RooNumIntConfig _convIntConfig
Configuration of numeric convolution integral ;.
RooAbsReal & model() const
RooAbsReal & pdf() const
TH2 * _callHist
! Histogram recording number of calls per convolution integral calculation
bool _useWindow
Switch to activate window convolution.
void setConvolutionWindow(RooAbsReal &centerParam, RooAbsReal &widthParam, double widthScaleFactor=1)
Restrict convolution integral to finite range [ x - C - S*W, x - C + S*W ] where x is current value o...
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
RooCategory & method1D()
RooCategory & method1DOpen()
std::unique_ptr< RooAbsIntegrator > createIntegrator(RooAbsFunc &func, const RooNumIntConfig &config, int ndim=0, bool isBinned=false) const
Construct a numeric integrator instance that operates on function 'func' and is configured with 'conf...
static RooNumIntFactory & instance()
Static method returning reference to singleton instance of factory.
static constexpr double infinity()
Return internal infinity representation.
Definition RooNumber.h:25
Variable that can be changed from the outside.
Definition RooRealVar.h:37
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
Int_t Fill(Double_t) override
Invalid Fill method.
Definition TH2.cxx:393
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Double_t x[n]
Definition legend1.C:17