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
76
77
78
79
80////////////////////////////////////////////////////////////////////////////////
81
83
84
85
86
87//_____________________________________________________________________________R
88RooNumConvPdf::RooNumConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& inPdf, RooAbsPdf& resmodel) :
89 RooAbsPdf(name,title),
90 _origVar("!origVar","Original Convolution variable",this,convVar),
91 _origPdf("!origPdf","Original Input PDF",this,inPdf),
92 _origModel("!origModel","Original Resolution model",this,resmodel)
93{
94 // Constructor of convolution operator PDF
95 //
96 // convVar : convolution variable (on which both pdf and resmodel should depend)
97 // pdf : input 'physics' pdf
98 // resmodel : input 'resolution' pdf
99 //
100 // output is pdf(x) (X) resmodel(x) = Int [ pdf(x') resmodel (x-x') ] dx'
101 //
102}
103
104
105
106////////////////////////////////////////////////////////////////////////////////
107/// Copy constructor
108
111 _origVar("!origVar",this,other._origVar),
112 _origPdf("!origPdf",this,other._origPdf),
113 _origModel("!origModel",this,other._origModel)
114{
115 // Make temporary clone of original convolution to preserve configuration information
116 // This information will be propagated to a newly create convolution in a subsequent
117 // call to initialize()
118 if (other._conv) {
119 _conv = std::make_unique<RooNumConvolution>(*other._conv,Form("%s_CONV",name?name:GetName())) ;
120 }
121}
122
123
124
125////////////////////////////////////////////////////////////////////////////////
126/// Destructor
127
129
130
131
132////////////////////////////////////////////////////////////////////////////////
133/// Calculate and return value of p.d.f
134
136{
137 if (!_init) initialize() ;
138
139 return _conv->evaluate() ;
140}
141
142
143
144////////////////////////////////////////////////////////////////////////////////
145/// One-time initialization of object
146
148{
149 // Save pointer to any prototype convolution object (only present if this object is made through
150 // a copy constructor)
152
153 // Optionally pass along configuration data from prototype object
154 _conv = std::make_unique<RooNumConvolution>(Form("%s_CONV",GetName()),GetTitle(),var(),pdf(),model(),protoConv) ;
155
156 _init = true ;
157}
158
159
160
161////////////////////////////////////////////////////////////////////////////////
162/// Return appropriate generator context for this convolved p.d.f. If both pdf and resolution
163/// model support internal generation return and optimization convolution generation context
164/// that uses a smearing algorithm. Otherwise return a standard accept/reject sampling
165/// context on the convoluted shape.
166
168 const RooArgSet* auxProto, bool verbose) const
169{
170 if (!_init) initialize() ;
171
172 // Check if physics PDF and resolution model can both directly generate the convolution variable
174 _conv->model().getObservables(&vars, modelDep);
175 modelDep.remove(_conv->var(),true,true) ;
176 Int_t numAddDep = modelDep.size() ;
177
178 RooArgSet dummy ;
179 bool pdfCanDir = ((static_cast<RooAbsPdf&>(_conv->pdf())).getGenerator(_conv->var(),dummy) != 0 && \
180 (static_cast<RooAbsPdf&>(_conv->pdf())).isDirectGenSafe(_conv->var())) ;
181 bool resCanDir = ((static_cast<RooAbsPdf&>(_conv->model())).getGenerator(_conv->var(),dummy) !=0 &&
182 (static_cast<RooAbsPdf&>(_conv->model())).isDirectGenSafe(_conv->var())) ;
183
184 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
185 // Any resolution model with more dependents than the convolution variable
186 // or pdf or resmodel do not support direct generation
187 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
188 }
189
190 // Any other resolution model: use specialized generator context
191 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
192}
193
194
195
196////////////////////////////////////////////////////////////////////////////////
197/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
198/// product operator construction
199
200void RooNumConvPdf::printMetaArgs(ostream& os) const
201{
202 os << _origPdf.arg().GetName() << "(" << _origVar.arg().GetName() << ") (*) " << _origModel.arg().GetName() << "(" << _origVar.arg().GetName() << ") " ;
203}
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
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:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50