Logo ROOT   6.10/09
Reference Guide
rf111_numintconfig.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'BASIC FUNCTIONALITY' RooFit tutorial macro #111
4 //
5 // Configuration and customization of how numeric (partial) integrals
6 // are executed
7 //
8 //
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooDataSet.h"
19 #include "RooGaussian.h"
20 #include "TCanvas.h"
21 #include "RooPlot.h"
22 #include "RooNumIntConfig.h"
23 #include "RooLandau.h"
24 #include "RooArgSet.h"
25 #include <iomanip>
26 using namespace RooFit ;
27 
28 class TestBasic111 : public RooFitTestUnit
29 {
30 public:
31  TestBasic111(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Numeric integration configuration",refFile,writeRef,verbose) {} ;
32  Bool_t testCode() {
33 
34  // 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
35  // ----------------------------------------------------------------------------
36 
37  // Example: Change global precision for 1D integrals from 1e-7 to 1e-6
38  //
39  // The relative epsilon (change as fraction of current best integral estimate) and
40  // absolute epsilon (absolute change w.r.t last best integral estimate) can be specified
41  // separately. For most p.d.f integrals the relative change criterium is the most important,
42  // however for certain non-p.d.f functions that integrate out to zero a separate absolute
43  // change criterium is necessary to declare convergence of the integral
44  //
45  // NB: This change is for illustration only. In general the precision should be at least 1e-7
46  // for normalization integrals for MINUIT to succeed.
47  //
50 
51 
52  // 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
53  // ------------------------------------------------------------------
54 
55  // Construct p.d.f without support for analytical integrator for demonstration purposes
56  RooRealVar x("x","x",-10,10) ;
57  RooLandau landau("landau","landau",x,RooConst(0),RooConst(0.1)) ;
58 
59 
60  // Calculate integral over landau with default choice of numeric integrator
61  RooAbsReal* intLandau = landau.createIntegral(x) ;
62  Double_t val = intLandau->getVal() ;
63  regValue(val,"rf111_val1") ;
64 
65 
66  // 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
67  // -----------------------------------------------------------
68 
69 
70  // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
71  // for closed 1D integrals
73  customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
74 
75 
76  // Calculate integral over landau with custom integral specification
77  RooAbsReal* intLandau2 = landau.createIntegral(x,NumIntConfig(customConfig)) ;
78  Double_t val2 = intLandau2->getVal() ;
79  regValue(val2,"rf111_val2") ;
80 
81 
82  // 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
83  // -------------------------------------------------------------------------------------
84 
85 
86  // Another possibility: associate custom numeric integration configuration as default for object 'landau'
87  landau.setIntegratorConfig(customConfig) ;
88 
89 
90  // Calculate integral over landau custom numeric integrator specified as object default
91  RooAbsReal* intLandau3 = landau.createIntegral(x) ;
92  Double_t val3 = intLandau3->getVal() ;
93  regValue(val3,"rf111_val3") ;
94 
95 
96  delete intLandau ;
97  delete intLandau2 ;
98  delete intLandau3 ;
99 
100  return kTRUE ;
101  }
102 } ;
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Landau Distribution p.d.f.
Definition: RooLandau.h:24
Double_t x[n]
Definition: legend1.C:17
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:501
bool verbose
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
RooConstVar & RooConst(Double_t val)
void setEpsAbs(Double_t newEpsAbs)
Set absolute convergence criteria (convergence if abs(Err)<newEpsAbs)
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooCmdArg NumIntConfig(const RooNumIntConfig &cfg)
void setEpsRel(Double_t newEpsRel)
Set relative convergence criteria (convergence if abs(Err)/abs(Int)<newEpsRel)