Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooNumConvPdf.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 RooNumConvPdf.cxx
19\class RooNumConvPdf
20\ingroup Roofitcore
21
22Numeric 1-dimensional convolution operator PDF. This class can convolve any PDF
23with any other PDF using a straightforward numeric calculation of the
24convolution integral
25This class should be used as last resort as numeric convolution calculated
26this way is computationally intensive and prone to stability fitting problems.
27<b>The preferred way to compute numeric convolutions is RooFFTConvPdf</b>,
28which calculates convolutions using Fourier Transforms (requires external free
29FFTW3 package)
30RooNumConvPdf implements reasonable defaults that should convolve most
31functions reasonably well, but results strongly depend on the shape of your
32input PDFS so always check your result.
33The default integration engine for the numeric convolution is the
34adaptive Gauss-Kronrod method, which empirically seems the most robust
35for this task. You can override the convolution integration settings via
36the RooNumIntConfig object reference returned by the convIntConfig() member
37function
38By default the numeric convolution is integrated from -infinity to
39+infinity through a <pre>x -> 1/x</pre> coordinate transformation of the
40tails. For convolution with a very small bandwidth it may be
41advantageous (for both CPU consumption and stability) if the
42integration domain is limited to a finite range. The function
43setConvolutionWindow(mean,width,scale) allows to set a sliding
44window around the x value to be calculated taking a RooAbsReal
45expression for an offset and a width to be taken around the x
46value. These input expression can be RooFormulaVars or other
47function objects although the 3d 'scale' argument 'scale'
48multiplies the width RooAbsReal expression given in the 2nd
49argument, allowing for an appropriate window definition for most
50cases without need for a RooFormulaVar object: e.g. a Gaussian
51resolution PDF do setConvolutionWindow(gaussMean,gaussSigma,5)
52Note that for a 'wide' Gaussian the -inf to +inf integration
53may converge more quickly than that over a finite range!
54The default numeric precision is 1e-7, i.e. the global default for
55numeric integration but you should experiment with this value to
56see if it is sufficient for example by studying the number of function
57calls that MINUIT needs to fit your function as function of the
58convolution precision.
59**/
60
61#include "Riostream.h"
62#include "RooNumConvPdf.h"
63#include "RooArgList.h"
64#include "RooRealVar.h"
65#include "RooFormulaVar.h"
66#include "RooCustomizer.h"
68#include "RooNumIntFactory.h"
69#include "RooGenContext.h"
70#include "RooConvGenContext.h"
71
72
73
74using std::ostream;
75
77
78
79
80
81////////////////////////////////////////////////////////////////////////////////
82
84
85
86
87
88//_____________________________________________________________________________R
89RooNumConvPdf::RooNumConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& inPdf, RooAbsPdf& resmodel) :
90 RooAbsPdf(name,title),
91 _origVar("!origVar","Original Convolution variable",this,convVar),
92 _origPdf("!origPdf","Original Input PDF",this,inPdf),
93 _origModel("!origModel","Original Resolution model",this,resmodel)
94{
95 // Constructor of convolution operator PDF
96 //
97 // convVar : convolution variable (on which both pdf and resmodel should depend)
98 // pdf : input 'physics' pdf
99 // resmodel : input 'resolution' pdf
100 //
101 // output is pdf(x) (X) resmodel(x) = Int [ pdf(x') resmodel (x-x') ] dx'
102 //
103}
104
105
106
107////////////////////////////////////////////////////////////////////////////////
108/// Copy constructor
109
111 RooAbsPdf(other,name),
112 _origVar("!origVar",this,other._origVar),
113 _origPdf("!origPdf",this,other._origPdf),
114 _origModel("!origModel",this,other._origModel)
115{
116 // Make temporary clone of original convolution to preserve configuration information
117 // This information will be propagated to a newly create convolution in a subsequent
118 // call to initialize()
119 if (other._conv) {
120 _conv = std::make_unique<RooNumConvolution>(*other._conv,Form("%s_CONV",name?name:GetName())) ;
121 }
122}
123
124
125
126////////////////////////////////////////////////////////////////////////////////
127/// Destructor
128
130
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Calculate and return value of p.d.f
135
137{
138 if (!_init) initialize() ;
139
140 return _conv->evaluate() ;
141}
142
143
144
145////////////////////////////////////////////////////////////////////////////////
146/// One-time initialization of object
147
149{
150 // Save pointer to any prototype convolution object (only present if this object is made through
151 // a copy constructor)
152 RooNumConvolution* protoConv = _conv.get();
153
154 // Optionally pass along configuration data from prototype object
155 _conv = std::make_unique<RooNumConvolution>(Form("%s_CONV",GetName()),GetTitle(),var(),pdf(),model(),protoConv) ;
156
157 _init = true ;
158}
159
160
161
162////////////////////////////////////////////////////////////////////////////////
163/// Return appropriate generator context for this convolved p.d.f. If both pdf and resolution
164/// model support internal generation return and optimization convolution generation context
165/// that uses a smearing algorithm. Otherwise return a standard accept/reject sampling
166/// context on the convoluted shape.
167
169 const RooArgSet* auxProto, bool verbose) const
170{
171 if (!_init) initialize() ;
172
173 // Check if physics PDF and resolution model can both directly generate the convolution variable
174 RooArgSet modelDep;
175 _conv->model().getObservables(&vars, modelDep);
176 modelDep.remove(_conv->var(),true,true) ;
177 Int_t numAddDep = modelDep.size() ;
178
179 RooArgSet dummy ;
180 bool pdfCanDir = ((static_cast<RooAbsPdf&>(_conv->pdf())).getGenerator(_conv->var(),dummy) != 0 && \
181 (static_cast<RooAbsPdf&>(_conv->pdf())).isDirectGenSafe(_conv->var())) ;
182 bool resCanDir = ((static_cast<RooAbsPdf&>(_conv->model())).getGenerator(_conv->var(),dummy) !=0 &&
183 (static_cast<RooAbsPdf&>(_conv->model())).isDirectGenSafe(_conv->var())) ;
184
185 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
186 // Any resolution model with more dependents than the convolution variable
187 // or pdf or resmodel do not support direct generation
188 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
189 }
190
191 // Any other resolution model: use specialized generator context
192 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
193}
194
195
196
197////////////////////////////////////////////////////////////////////////////////
198/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
199/// product operator construction
200
201void RooNumConvPdf::printMetaArgs(ostream& os) const
202{
203 os << _origPdf.arg().GetName() << "(" << _origVar.arg().GetName() << ") (*) " << _origModel.arg().GetName() << "(" << _origVar.arg().GetName() << ") " ;
204}
#define ClassImp(name)
Definition Rtypes.h:382
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
Storage_t::size_type size() const
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold unbinned data.
Definition RooDataSet.h:34
Implements a universal generator context for all RooAbsPdf classes that do not have or need a special...
Numeric 1-dimensional convolution operator PDF.
friend class RooConvGenContext
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Return appropriate generator context for this convolved p.d.f.
RooRealProxy _origPdf
Original input PDF.
~RooNumConvPdf() override
Destructor.
std::unique_ptr< RooNumConvolution > _conv
! Actual convolution calculation
RooRealProxy _origModel
Original resolution model.
RooRealVar & var() const
double evaluate() const override
Calculate and return value of p.d.f.
RooRealProxy _origVar
Original convolution variable.
RooAbsReal & pdf() const
RooAbsReal & model() const
bool _init
! do not persist
void printMetaArgs(std::ostream &os) const override
Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the p...
void initialize() const
One-time initialization of object.
Numeric 1-dimensional convolution operator PDF.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const T & arg() const
Return reference to object held in proxy.
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