Logo ROOT   6.14/05
Reference Guide
rf901_numintconfig.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'NUMERIC ALGORITHM TUNING' RooFit tutorial macro #901
5 ///
6 /// Configuration and customization of how numeric (partial) integrals
7 /// are executed
8 ///
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooConstVar.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 #include "RooNumIntConfig.h"
22 #include "RooLandau.h"
23 #include "RooArgSet.h"
24 #include <iomanip>
25 using namespace RooFit ;
26 
27 
28 void rf901_numintconfig()
29 {
30 
31  // 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
32  // ----------------------------------------------------------------------------
33 
34  // Print current global default configuration for numeric integration strategies
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  // Activate debug-level messages for topic integration to be able to follow actions below
62 
63 
64  // Calculate integral over landau with default choice of numeric integrator
65  RooAbsReal* intLandau = landau.createIntegral(x) ;
66  Double_t val = intLandau->getVal() ;
67  cout << " [1] int_dx landau(x) = " << setprecision(15) << val << endl ;
68 
69 
70 
71  // 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
72  // -----------------------------------------------------------
73 
74 
75  // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
76  // for closed 1D integrals
78  customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
79 
80 
81  // Calculate integral over landau with custom integral specification
82  RooAbsReal* intLandau2 = landau.createIntegral(x,NumIntConfig(customConfig)) ;
83  Double_t val2 = intLandau2->getVal() ;
84  cout << " [2] int_dx landau(x) = " << val2 << endl ;
85 
86 
87 
88  // 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
89  // -------------------------------------------------------------------------------------
90 
91 
92  // Another possibility: associate custom numeric integration configuration as default for object 'landau'
93  landau.setIntegratorConfig(customConfig) ;
94 
95 
96  // Calculate integral over landau custom numeric integrator specified as object default
97  RooAbsReal* intLandau3 = landau.createIntegral(x) ;
98  Double_t val3 = intLandau3->getVal() ;
99  cout << " [3] int_dx landau(x) = " << val3 << endl ;
100 
101 
102  // Another possibility: Change global default for 1D numeric integration strategy on finite domains
103  RooAbsReal::defaultIntegratorConfig()->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
104 
105 
106 
107  // 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
108  // ---------------------------------------------------------------------------------------
109 
110  // Adjust maximum number of steps of RooIntegrator1D in the global default configuration
111  RooAbsReal::defaultIntegratorConfig()->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",30) ;
112 
113 
114  // Example of how to change the parameters of a numeric integrator
115  // (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
116  // and RooCategories holding parameters with a finite set of options)
117  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg",50) ;
118  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method","15Points") ;
119 
120 
121  // Example of how to print set of possible values for "method" category
122  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method")->Print("v") ;
123 
124 }
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
Definition: RooArgSet.cxx:547
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
static RooMsgService & instance()
Return reference to singleton instance.
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)
Set value by specifying the name of the desired state If printError is set, a message will be printed...
Landau Distribution p.d.f.
Definition: RooLandau.h:24
RooCategory & method1D()
Double_t x[n]
Definition: legend1.C:17
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
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
RooCmdArg Topic(Int_t topic)
Int_t addStream(RooFit::MsgLevel level, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg())
Add a message logging stream for message with given RooFit::MsgLevel or higher (i.e.
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)
RooCmdArg NumIntConfig(const RooNumIntConfig &cfg)
void setEpsRel(Double_t newEpsRel)
Set relative convergence criteria (convergence if abs(Err)/abs(Int)<newEpsRel)