RooNumConvPdf.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooNumConvPdf.cxx 28259 2009-04-16 16:21:16Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // BEGIN_HTML
00020 // Numeric 1-dimensional convolution operator PDF. This class can convolve any PDF
00021 // with any other PDF using a straightforward numeric calculation of the
00022 // convolution integral
00023 // <p>
00024 // This class should be used as last resort as numeric convolution calculated
00025 // this way is computationally intensive and prone to stability fitting problems. 
00026 // <b>The preferred way to compute numeric convolutions is RooFFTConvPdf</b>,
00027 // which calculates convolutions using Fourier Transforms (requires external free
00028 // FFTW3 package)
00029 // <p>
00030 // RooNumConvPdf implements reasonable defaults that should convolve most
00031 // functions reasonably well, but results strongly depend on the shape of your
00032 // input PDFS so always check your result.
00033 // <p>
00034 // The default integration engine for the numeric convolution is the
00035 // adaptive Gauss-Kronrod method, which empirically seems the most robust
00036 // for this task. You can override the convolution integration settings via
00037 // the RooNumIntConfig object reference returned by the convIntConfig() member
00038 // function
00039 // <p>
00040 // By default the numeric convolution is integrated from -infinity to
00041 // +infinity through a <pre>x -> 1/x</pre> coordinate transformation of the
00042 // tails. For convolution with a very small bandwidth it may be
00043 // advantageous (for both CPU consumption and stability) if the
00044 // integration domain is limited to a finite range. The function
00045 // setConvolutionWindow(mean,width,scale) allows to set a sliding
00046 // window around the x value to be calculated taking a RooAbsReal
00047 // expression for an offset and a width to be taken around the x
00048 // value. These input expression can be RooFormulaVars or other
00049 // function objects although the 3d 'scale' argument 'scale'
00050 // multiplies the width RooAbsReal expression given in the 2nd
00051 // argument, allowing for an appropriate window definition for most
00052 // cases without need for a RooFormulaVar object: e.g. a Gaussian
00053 // resolution PDF do setConvolutionWindow(gaussMean,gaussSigma,5)
00054 // Note that for a 'wide' Gaussian the -inf to +inf integration
00055 // may converge more quickly than that over a finite range!
00056 // <p>
00057 // The default numeric precision is 1e-7, i.e. the global default for
00058 // numeric integration but you should experiment with this value to
00059 // see if it is sufficient for example by studying the number of function
00060 // calls that MINUIT needs to fit your function as function of the
00061 // convolution precision. 
00062 // END_HTML
00063 //
00064 
00065 #include "RooFit.h"
00066 
00067 #include "Riostream.h"
00068 #include "Riostream.h"
00069 #include "TH2F.h"
00070 #include "RooNumConvPdf.h"
00071 #include "RooArgList.h"
00072 #include "RooRealVar.h"
00073 #include "RooFormulaVar.h"
00074 #include "RooCustomizer.h"
00075 #include "RooConvIntegrandBinding.h"
00076 #include "RooNumIntFactory.h"
00077 #include "RooGenContext.h"
00078 #include "RooConvGenContext.h"
00079 
00080 
00081 
00082 ClassImp(RooNumConvPdf)
00083 ;
00084 
00085 
00086 
00087 //_____________________________________________________________________________
00088 RooNumConvPdf::RooNumConvPdf() :
00089   _init(kFALSE),
00090   _conv(0)
00091 {  
00092 }
00093 
00094 
00095 
00096 
00097 //_____________________________________________________________________________R
00098 RooNumConvPdf::RooNumConvPdf(const char *name, const char *title, RooRealVar& convVar, RooAbsPdf& inPdf, RooAbsPdf& resmodel) : 
00099   RooAbsPdf(name,title), 
00100   _init(kFALSE),
00101   _conv(0),
00102   _origVar("!origVar","Original Convolution variable",this,convVar),
00103   _origPdf("!origPdf","Original Input PDF",this,inPdf),
00104   _origModel("!origModel","Original Resolution model",this,resmodel)
00105 {
00106   // Constructor of convolution operator PDF
00107   // 
00108   // convVar  :  convolution variable (on which both pdf and resmodel should depend)
00109   // pdf      :  input 'physics' pdf
00110   // resmodel :  input 'resultion' pdf
00111   //
00112   // output is pdf(x) (X) resmodel(x) = Int [ pdf(x') resmodel (x-x') ] dx'
00113   //
00114 }
00115 
00116 
00117 
00118 //_____________________________________________________________________________
00119 RooNumConvPdf::RooNumConvPdf(const RooNumConvPdf& other, const char* name) :
00120   RooAbsPdf(other,name), 
00121   _init(kFALSE),
00122   _origVar("!origVar",this,other._origVar),
00123   _origPdf("!origPdf",this,other._origPdf),
00124   _origModel("!origModel",this,other._origModel)
00125 {
00126   // Copy constructor
00127 
00128   // Make temporary clone of original convolution to preserve configuration information
00129   // This information will be propagated to a newly create convolution in a subsequent
00130   // call to initialize() 
00131   if (other._conv) {
00132     _conv = new RooNumConvolution(*other._conv,Form("%s_CONV",name?name:GetName())) ;
00133   } else {
00134     _conv = 0 ;
00135   }
00136 }
00137 
00138 
00139 
00140 //_____________________________________________________________________________
00141 RooNumConvPdf::~RooNumConvPdf() 
00142 {
00143   // Destructor
00144   if (_init) {
00145     delete _conv ;
00146   }
00147 }
00148 
00149 
00150 
00151 //_____________________________________________________________________________
00152 Double_t RooNumConvPdf::evaluate() const
00153 {
00154   // Calculate and return value of p.d.f
00155 
00156   if (!_init) initialize() ;
00157 
00158   return _conv->evaluate() ;
00159 }
00160 
00161 
00162 
00163 //_____________________________________________________________________________
00164 void RooNumConvPdf::initialize() const
00165 {
00166   // One-time initialization of object
00167 
00168   // Save pointer to any prototype convolution object (only present if this object is made through
00169   // a copy constructor) 
00170   RooNumConvolution* protoConv = _conv ;
00171 
00172   // Optionally pass along configuration data from prototype object
00173   _conv = new RooNumConvolution(Form("%s_CONV",GetName()),GetTitle(),var(),pdf(),model(),protoConv) ;
00174 
00175   // Delete prototype object now
00176   if (protoConv) {
00177     delete protoConv ;  
00178   }
00179 
00180   _init = kTRUE ;
00181 }
00182 
00183 
00184 
00185 //_____________________________________________________________________________
00186 RooAbsGenContext* RooNumConvPdf::genContext(const RooArgSet &vars, const RooDataSet *prototype, 
00187                                             const RooArgSet* auxProto, Bool_t verbose) const 
00188 {
00189   // Return appropriate generator context for this convolved p.d.f. If both pdf and resolution
00190   // model support internal generation return and optimization convolution generation context
00191   // that uses a smearing algorithm. Otherwise return a standard accept/reject sampling
00192   // context on the convoluted shape.
00193 
00194   if (!_init) initialize() ;
00195 
00196   // Check if physics PDF and resolution model can both directly generate the convolution variable
00197   RooArgSet* modelDep = _conv->model().getObservables(&vars) ;
00198   modelDep->remove(_conv->var(),kTRUE,kTRUE) ;
00199   Int_t numAddDep = modelDep->getSize() ;
00200   delete modelDep ;
00201 
00202   RooArgSet dummy ;
00203   Bool_t pdfCanDir = (((RooAbsPdf&)_conv->pdf()).getGenerator(_conv->var(),dummy) != 0 && \
00204                       ((RooAbsPdf&)_conv->pdf()).isDirectGenSafe(_conv->var())) ;
00205   Bool_t resCanDir = (((RooAbsPdf&)_conv->model()).getGenerator(_conv->var(),dummy) !=0  && 
00206                       ((RooAbsPdf&)_conv->model()).isDirectGenSafe(_conv->var())) ;
00207 
00208   if (numAddDep>0 || !pdfCanDir || !resCanDir) {
00209     // Any resolution model with more dependents than the convolution variable
00210     // or pdf or resmodel do not support direct generation
00211     return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
00212   } 
00213   
00214   // Any other resolution model: use specialized generator context
00215   return new RooConvGenContext(*this,vars,prototype,auxProto,verbose) ;
00216 }
00217 
00218 
00219 
00220 //_____________________________________________________________________________
00221 void RooNumConvPdf::printMetaArgs(ostream& os) const 
00222 {
00223   // Customized printing of arguments of a RooNumConvPdf to more intuitively reflect the contents of the
00224   // product operator construction
00225 
00226   os << _origPdf.arg().GetName() << "(" << _origVar.arg().GetName() << ") (*) " << _origModel.arg().GetName() << "(" << _origVar.arg().GetName() << ") " ;
00227 }

Generated on Tue Jul 5 15:06:53 2011 for ROOT_528-00b_version by  doxygen 1.5.1