Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf901_numintconfig.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// Numeric algorithm tuning: configuration and customization of how numeric (partial) integrals are executed
5 ///
6 /// \macro_output
7 /// \macro_code
8 /// \author 07/2008 - Wouter Verkerke
9 
10 #include "RooRealVar.h"
11 #include "RooDataSet.h"
12 #include "RooGaussian.h"
13 #include "RooConstVar.h"
14 #include "TCanvas.h"
15 #include "TAxis.h"
16 #include "RooPlot.h"
17 #include "RooNumIntConfig.h"
18 #include "RooLandau.h"
19 #include "RooArgSet.h"
20 #include <iomanip>
21 using namespace RooFit;
22 
23 void rf901_numintconfig()
24 {
25 
26  // 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
27  // ----------------------------------------------------------------------------
28 
29  // Print current global default configuration for numeric integration strategies
30  RooAbsReal::defaultIntegratorConfig()->Print("v");
31 
32  // Example: Change global precision for 1D integrals from 1e-7 to 1e-6
33  //
34  // The relative epsilon (change as fraction of current best integral estimate) and
35  // absolute epsilon (absolute change w.r.t last best integral estimate) can be specified
36  // separately. For most p.d.f integrals the relative change criterium is the most important,
37  // however for certain non-p.d.f functions that integrate out to zero a separate absolute
38  // change criterium is necessary to declare convergence of the integral
39  //
40  // NB: This change is for illustration only. In general the precision should be at least 1e-7
41  // for normalization integrals for MINUIT to succeed.
42  //
43  RooAbsReal::defaultIntegratorConfig()->setEpsAbs(1e-6);
44  RooAbsReal::defaultIntegratorConfig()->setEpsRel(1e-6);
45 
46  // 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
47  // ------------------------------------------------------------------
48 
49  // Construct p.d.f without support for analytical integrator for demonstration purposes
50  RooRealVar x("x", "x", -10, 10);
51  RooLandau landau("landau", "landau", x, RooConst(0), RooConst(0.1));
52 
53  // Activate debug-level messages for topic integration to be able to follow actions below
54  RooMsgService::instance().addStream(DEBUG, Topic(Integration));
55 
56  // Calculate integral over landau with default choice of numeric integrator
57  RooAbsReal *intLandau = landau.createIntegral(x);
58  Double_t val = intLandau->getVal();
59  cout << " [1] int_dx landau(x) = " << setprecision(15) << val << endl;
60 
61  // 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
62  // -----------------------------------------------------------
63 
64  // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
65  // for closed 1D integrals
66  RooNumIntConfig customConfig(*RooAbsReal::defaultIntegratorConfig());
67 #ifdef R__HAS_MATHMORE
68  customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D");
69 #else
70  Warning("rf901_numintconfig","ROOT is built without Mathmore (GSL) support. Cannot use RooAdaptiveGaussKronrodIntegrator1D");
71 #endif
72 
73  // Calculate integral over landau with custom integral specification
74  RooAbsReal *intLandau2 = landau.createIntegral(x, NumIntConfig(customConfig));
75  Double_t val2 = intLandau2->getVal();
76  cout << " [2] int_dx landau(x) = " << val2 << endl;
77 
78  // 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
79  // -------------------------------------------------------------------------------------
80 
81  // Another possibility: associate custom numeric integration configuration as default for object 'landau'
82  landau.setIntegratorConfig(customConfig);
83 
84  // Calculate integral over landau custom numeric integrator specified as object default
85  RooAbsReal *intLandau3 = landau.createIntegral(x);
86  Double_t val3 = intLandau3->getVal();
87  cout << " [3] int_dx landau(x) = " << val3 << endl;
88 
89  // Another possibility: Change global default for 1D numeric integration strategy on finite domains
90 #ifdef R__HAS_MATHMORE
91  RooAbsReal::defaultIntegratorConfig()->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D");
92 
93  // 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
94  // ---------------------------------------------------------------------------------------
95 
96  // Adjust maximum number of steps of RooIntegrator1D in the global default configuration
97  RooAbsReal::defaultIntegratorConfig()->getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 30);
98 
99  // Example of how to change the parameters of a numeric integrator
100  // (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
101  // and RooCategories holding parameters with a finite set of options)
102  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg", 50);
103  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method", "15Points");
104 
105  // Example of how to print set of possible values for "method" category
106  customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method")->Print("v");
107 #endif
108 
109 }