rf901_numintconfig.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'NUMERIC ALGORITHM TUNING' RooFit tutorial macro #901 
00004 // 
00005 // Configuration and customization of how numeric (partial) integrals
00006 // are executed
00007 //
00008 //
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 // 
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooConstVar.h"
00021 #include "TCanvas.h"
00022 #include "TAxis.h"
00023 #include "RooPlot.h"
00024 #include "RooNumIntConfig.h"
00025 #include "RooLandau.h"
00026 #include "RooArgSet.h"
00027 #include <iomanip>
00028 using namespace RooFit ;
00029 
00030 
00031 void rf901_numintconfig()
00032 {
00033 
00034   // A d j u s t   g l o b a l   1 D   i n t e g r a t i o n   p r e c i s i o n 
00035   // ----------------------------------------------------------------------------
00036 
00037   // Print current global default configuration for numeric integration strategies
00038   RooAbsReal::defaultIntegratorConfig()->Print("v") ;
00039 
00040   // Example: Change global precision for 1D integrals from 1e-7 to 1e-6
00041   //
00042   // The relative epsilon (change as fraction of current best integral estimate) and
00043   // absolute epsilon (absolute change w.r.t last best integral estimate) can be specified
00044   // separately. For most p.d.f integrals the relative change criterium is the most important,
00045   // however for certain non-p.d.f functions that integrate out to zero a separate absolute
00046   // change criterium is necessary to declare convergence of the integral
00047   //
00048   // NB: This change is for illustration only. In general the precision should be at least 1e-7 
00049   // for normalization integrals for MINUIT to succeed.
00050   //
00051   RooAbsReal::defaultIntegratorConfig()->setEpsAbs(1e-6) ;
00052   RooAbsReal::defaultIntegratorConfig()->setEpsRel(1e-6) ;
00053 
00054 
00055   // N u m e r i c   i n t e g r a t i o n   o f   l a n d a u   p d f 
00056   // ------------------------------------------------------------------
00057   
00058   // Construct p.d.f without support for analytical integrator for demonstration purposes
00059   RooRealVar x("x","x",-10,10) ;
00060   RooLandau landau("landau","landau",x,RooConst(0),RooConst(0.1)) ;
00061   
00062 
00063   // Activate debug-level messages for topic integration to be able to follow actions below
00064   RooMsgService::instance().addStream(DEBUG,Topic(Integration)) ;
00065 
00066 
00067   // Calculate integral over landau with default choice of numeric integrator
00068   RooAbsReal* intLandau = landau.createIntegral(x) ;
00069   Double_t val = intLandau->getVal() ;
00070   cout << " [1] int_dx landau(x) = " << setprecision(15) << val << endl ;
00071 
00072 
00073 
00074   // S a m e   w i t h   c u s t o m   c o n f i g u r a t i o n
00075   // -----------------------------------------------------------
00076   
00077 
00078   // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
00079   // for closed 1D integrals
00080   RooNumIntConfig customConfig(*RooAbsReal::defaultIntegratorConfig()) ;
00081   customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
00082 
00083 
00084   // Calculate integral over landau with custom integral specification
00085   RooAbsReal* intLandau2 = landau.createIntegral(x,NumIntConfig(customConfig)) ;
00086   Double_t val2 = intLandau2->getVal() ;
00087   cout << " [2] int_dx landau(x) = " << val2 << endl ;
00088 
00089 
00090 
00091   // A d j u s t i n g   d e f a u l t   c o n f i g   f o r   a   s p e c i f i c   p d f 
00092   // -------------------------------------------------------------------------------------
00093   
00094 
00095   // Another possibility: associate custom numeric integration configuration as default for object 'landau'
00096   landau.setIntegratorConfig(customConfig) ;
00097 
00098 
00099   // Calculate integral over landau custom numeric integrator specified as object default
00100   RooAbsReal* intLandau3 = landau.createIntegral(x) ;
00101   Double_t val3 = intLandau3->getVal() ;
00102   cout << " [3] int_dx landau(x) = " << val3 << endl ;
00103  
00104 
00105   // Another possibility: Change global default for 1D numeric integration strategy on finite domains
00106   RooAbsReal::defaultIntegratorConfig()->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;  
00107 
00108 
00109 
00110   // A d j u s t i n g   p a r a m e t e r s   o f   a   s p e c i f i c   t e c h n i q u e 
00111   // ---------------------------------------------------------------------------------------
00112 
00113   // Adjust maximum number of steps of RooIntegrator1D in the global default configuration
00114   RooAbsReal::defaultIntegratorConfig()->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",30) ;
00115 
00116  
00117   // Example of how to change the parameters of a numeric integrator
00118   // (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
00119   //  and RooCategories holding parameters with a finite set of options)
00120   customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg",50) ;
00121   customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method","15Points") ;
00122 
00123 
00124   // Example of how to print set of possible values for "method" category
00125   customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method")->Print("v") ;
00126 
00127 }

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