00001 ////////////////////////////////////////////////////////////////////////// 00002 // 00003 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #111 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 "TCanvas.h" 00021 #include "RooPlot.h" 00022 #include "RooNumIntConfig.h" 00023 #include "RooLandau.h" 00024 #include "RooArgSet.h" 00025 #include <iomanip> 00026 using namespace RooFit ; 00027 00028 class TestBasic111 : public RooFitTestUnit 00029 { 00030 public: 00031 TestBasic111(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Numeric integration configuration",refFile,writeRef,verbose) {} ; 00032 Bool_t testCode() { 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 // Example: Change global precision for 1D integrals from 1e-7 to 1e-6 00038 // 00039 // The relative epsilon (change as fraction of current best integral estimate) and 00040 // absolute epsilon (absolute change w.r.t last best integral estimate) can be specified 00041 // separately. For most p.d.f integrals the relative change criterium is the most important, 00042 // however for certain non-p.d.f functions that integrate out to zero a separate absolute 00043 // change criterium is necessary to declare convergence of the integral 00044 // 00045 // NB: This change is for illustration only. In general the precision should be at least 1e-7 00046 // for normalization integrals for MINUIT to succeed. 00047 // 00048 RooAbsReal::defaultIntegratorConfig()->setEpsAbs(1e-6) ; 00049 RooAbsReal::defaultIntegratorConfig()->setEpsRel(1e-6) ; 00050 00051 00052 // 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 00053 // ------------------------------------------------------------------ 00054 00055 // Construct p.d.f without support for analytical integrator for demonstration purposes 00056 RooRealVar x("x","x",-10,10) ; 00057 RooLandau landau("landau","landau",x,RooConst(0),RooConst(0.1)) ; 00058 00059 00060 // Calculate integral over landau with default choice of numeric integrator 00061 RooAbsReal* intLandau = landau.createIntegral(x) ; 00062 Double_t val = intLandau->getVal() ; 00063 regValue(val,"rf111_val1") ; 00064 00065 00066 // 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 00067 // ----------------------------------------------------------- 00068 00069 00070 // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique 00071 // for closed 1D integrals 00072 RooNumIntConfig customConfig(*RooAbsReal::defaultIntegratorConfig()) ; 00073 customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ; 00074 00075 00076 // Calculate integral over landau with custom integral specification 00077 RooAbsReal* intLandau2 = landau.createIntegral(x,NumIntConfig(customConfig)) ; 00078 Double_t val2 = intLandau2->getVal() ; 00079 regValue(val2,"rf111_val2") ; 00080 00081 00082 // 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 00083 // ------------------------------------------------------------------------------------- 00084 00085 00086 // Another possibility: associate custom numeric integration configuration as default for object 'landau' 00087 landau.setIntegratorConfig(customConfig) ; 00088 00089 00090 // Calculate integral over landau custom numeric integrator specified as object default 00091 RooAbsReal* intLandau3 = landau.createIntegral(x) ; 00092 Double_t val3 = intLandau3->getVal() ; 00093 regValue(val3,"rf111_val3") ; 00094 00095 00096 delete intLandau ; 00097 delete intLandau2 ; 00098 delete intLandau3 ; 00099 00100 return kTRUE ; 00101 } 00102 } ;