Logo ROOT  
Reference Guide
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 "TH2F.h"
63#include "RooNumConvPdf.h"
64#include "RooArgList.h"
65#include "RooRealVar.h"
66#include "RooFormulaVar.h"
67#include "RooCustomizer.h"
69#include "RooNumIntFactory.h"
70#include "RooGenContext.h"
71#include "RooConvGenContext.h"
72
73
74
75using namespace std;
76
78
79
80
81
82////////////////////////////////////////////////////////////////////////////////
83
85 _init(kFALSE),
86 _conv(0)
87{
88}
89
90
91
92
93//_____________________________________________________________________________R
94RooNumConvPdf::RooNumConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& inPdf, RooAbsPdf& resmodel) :
95 RooAbsPdf(name,title),
96 _init(kFALSE),
97 _conv(0),
98 _origVar("!origVar","Original Convolution variable",this,convVar),
99 _origPdf("!origPdf","Original Input PDF",this,inPdf),
100 _origModel("!origModel","Original Resolution model",this,resmodel)
101{
102 // Constructor of convolution operator PDF
103 //
104 // convVar : convolution variable (on which both pdf and resmodel should depend)
105 // pdf : input 'physics' pdf
106 // resmodel : input 'resultion' pdf
107 //
108 // output is pdf(x) (X) resmodel(x) = Int [ pdf(x') resmodel (x-x') ] dx'
109 //
110}
111
112
113
114////////////////////////////////////////////////////////////////////////////////
115/// Copy constructor
116
118 RooAbsPdf(other,name),
119 _init(kFALSE),
120 _origVar("!origVar",this,other._origVar),
121 _origPdf("!origPdf",this,other._origPdf),
122 _origModel("!origModel",this,other._origModel)
123{
124 // Make temporary clone of original convolution to preserve configuration information
125 // This information will be propagated to a newly create convolution in a subsequent
126 // call to initialize()
127 if (other._conv) {
128 _conv = new RooNumConvolution(*other._conv,Form("%s_CONV",name?name:GetName())) ;
129 } else {
130 _conv = 0 ;
131 }
132}
133
134
135
136////////////////////////////////////////////////////////////////////////////////
137/// Destructor
138
140{
141 if (_init) {
142 delete _conv ;
143 }
144}
145
146
147
148////////////////////////////////////////////////////////////////////////////////
149/// Calculate and return value of p.d.f
150
152{
153 if (!_init) initialize() ;
154
155 return _conv->evaluate() ;
156}
157
158
159
160////////////////////////////////////////////////////////////////////////////////
161/// One-time initialization of object
162
164{
165 // Save pointer to any prototype convolution object (only present if this object is made through
166 // a copy constructor)
167 RooNumConvolution* protoConv = _conv ;
168
169 // Optionally pass along configuration data from prototype object
170 _conv = new RooNumConvolution(Form("%s_CONV",GetName()),GetTitle(),var(),pdf(),model(),protoConv) ;
171
172 // Delete prototype object now
173 if (protoConv) {
174 delete protoConv ;
175 }
176
177 _init = kTRUE ;
178}
179
180
181
182////////////////////////////////////////////////////////////////////////////////
183/// Return appropriate generator context for this convolved p.d.f. If both pdf and resolution
184/// model support internal generation return and optimization convolution generation context
185/// that uses a smearing algorithm. Otherwise return a standard accept/reject sampling
186/// context on the convoluted shape.
187
189 const RooArgSet* auxProto, Bool_t verbose) const
190{
191 if (!_init) initialize() ;
192
193 // Check if physics PDF and resolution model can both directly generate the convolution variable
194 RooArgSet* modelDep = _conv->model().getObservables(&vars) ;
195 modelDep->remove(_conv->var(),kTRUE,kTRUE) ;
196 Int_t numAddDep = modelDep->getSize() ;
197 delete modelDep ;
198
199 RooArgSet dummy ;
200 Bool_t pdfCanDir = (((RooAbsPdf&)_conv->pdf()).getGenerator(_conv->var(),dummy) != 0 && \
202 Bool_t resCanDir = (((RooAbsPdf&)_conv->model()).getGenerator(_conv->var(),dummy) !=0 &&
204
205 if (numAddDep>0 || !pdfCanDir || !resCanDir) {
206 // Any resolution model with more dependents than the convolution variable
207 // or pdf or resmodel do not support direct generation
208 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
209 }
210
211 // Any other resolution model: use specialized generator context
212 return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
213}
214
215
216
217////////////////////////////////////////////////////////////////////////////////
218/// Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
219/// product operator construction
220
221void RooNumConvPdf::printMetaArgs(ostream& os) const
222{
223 os << _origPdf.arg().GetName() << "(" << _origVar.arg().GetName() << ") (*) " << _origModel.arg().GetName() << "(" << _origVar.arg().GetName() << ") " ;
224}
const Bool_t kFALSE
Definition: RtypesCore.h:101
const Bool_t kTRUE
Definition: RtypesCore.h:100
#define ClassImp(name)
Definition: Rtypes.h:375
char name[80]
Definition: TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition: TString.cxx:2447
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition: RooAbsArg.h:318
Int_t getSize() const
Return the number of elements in the collection.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition: RooAbsPdf.cxx:2302
virtual Bool_t isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Definition: RooAbsPdf.cxx:2337
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:57
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:55
Class RooGenContext implement a universal generator context for all RooAbsPdf classes that do not hav...
Definition: RooGenContext.h:30
Numeric 1-dimensional convolution operator PDF.
Definition: RooNumConvPdf.h:26
friend class RooConvGenContext
Definition: RooNumConvPdf.h:76
RooRealProxy _origPdf
Original input PDF.
Definition: RooNumConvPdf.h:70
~RooNumConvPdf() override
Destructor.
RooNumConvolution * _conv
! Actual convolution calculation
Definition: RooNumConvPdf.h:67
Bool_t _init
! do not persist
Definition: RooNumConvPdf.h:65
RooRealProxy _origModel
Original resolution model.
Definition: RooNumConvPdf.h:71
RooRealVar & var() const
Definition: RooNumConvPdf.h:52
RooRealProxy _origVar
Original convolution variable.
Definition: RooNumConvPdf.h:69
RooAbsReal & pdf() const
Definition: RooNumConvPdf.h:53
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const override
Return appropriate generator context for this convolved p.d.f.
RooAbsReal & model() const
Definition: RooNumConvPdf.h:54
Double_t evaluate() const override
Calculate and return value of p.d.f.
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.
RooRealVar & var() const
RooAbsReal & model() const
RooAbsReal & pdf() const
Double_t evaluate() const override
Calculate convolution integral.
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:40
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