Logo ROOT  
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: 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>
21using namespace RooFit;
22
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
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 //
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
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
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}
#define DEBUG
Definition: Polynomial.cxx:41
#define e(i)
Definition: RSha256.hxx:103
double Double_t
Definition: RtypesCore.h:55
void Warning(const char *location, const char *msgfmt,...)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
static RooNumIntConfig * defaultIntegratorConfig()
Returns the default numeric integration configuration for all RooAbsReals.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
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:562
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:493
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
static RooMsgService & instance()
Return reference to singleton instance.
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....
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
void setEpsAbs(Double_t newEpsAbs)
Set absolute convergence criteria (convergence if abs(Err)<newEpsAbs)
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
void setEpsRel(Double_t newEpsRel)
Set relative convergence criteria (convergence if abs(Err)/abs(Int)<newEpsRel)
RooCategory & method1D()
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg Topic(Int_t topic)
RooCmdArg NumIntConfig(const RooNumIntConfig &cfg)
@ Integration
Definition: RooGlobalFunc.h:67
RooConstVar & RooConst(Double_t val)