Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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///
9/// \date July 2008
10/// \author Wouter Verkerke
11
12#include "RooRealVar.h"
13#include "RooDataSet.h"
14#include "RooGaussian.h"
15#include "RooConstVar.h"
16#include "TCanvas.h"
17#include "TAxis.h"
18#include "RooPlot.h"
19#include "RooNumIntConfig.h"
20#include "RooLandau.h"
21#include "RooArgSet.h"
22#include <iomanip>
23using namespace RooFit;
24
26{
27
28 // 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
29 // ----------------------------------------------------------------------------
30
31 // Print current global default configuration for numeric integration strategies
33
34 // Example: Change global precision for 1D integrals from 1e-7 to 1e-6
35 //
36 // The relative epsilon (change as fraction of current best integral estimate) and
37 // absolute epsilon (absolute change w.r.t last best integral estimate) can be specified
38 // separately. For most pdf integrals the relative change criterium is the most important,
39 // however for certain non-pdf functions that integrate out to zero a separate absolute
40 // change criterium is necessary to declare convergence of the integral
41 //
42 // NB: This change is for illustration only. In general the precision should be at least 1e-7
43 // for normalization integrals for MINUIT to succeed.
44 //
47
48 // 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
49 // ------------------------------------------------------------------
50
51 // Construct pdf without support for analytical integrator for demonstration purposes
52 RooRealVar x("x", "x", -10, 10);
53 RooLandau landau("landau", "landau", x, RooConst(0), RooConst(0.1));
54
55 // Activate debug-level messages for topic integration to be able to follow actions below
57
58 // Calculate integral over landau with default choice of numeric integrator
59 RooAbsReal *intLandau = landau.createIntegral(x);
60 Double_t val = intLandau->getVal();
61 cout << " [1] int_dx landau(x) = " << setprecision(15) << val << endl;
62
63 // 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
64 // -----------------------------------------------------------
65
66 // Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
67 // for closed 1D integrals
69#ifdef R__HAS_MATHMORE
70 customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D");
71#else
72 Warning("rf901_numintconfig","ROOT is built without Mathmore (GSL) support. Cannot use RooAdaptiveGaussKronrodIntegrator1D");
73#endif
74
75 // Calculate integral over landau with custom integral specification
76 RooAbsReal *intLandau2 = landau.createIntegral(x, NumIntConfig(customConfig));
77 Double_t val2 = intLandau2->getVal();
78 cout << " [2] int_dx landau(x) = " << val2 << endl;
79
80 // 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
81 // -------------------------------------------------------------------------------------
82
83 // Another possibility: associate custom numeric integration configuration as default for object 'landau'
84 landau.setIntegratorConfig(customConfig);
85
86 // Calculate integral over landau custom numeric integrator specified as object default
87 RooAbsReal *intLandau3 = landau.createIntegral(x);
88 Double_t val3 = intLandau3->getVal();
89 cout << " [3] int_dx landau(x) = " << val3 << endl;
90
91 // Another possibility: Change global default for 1D numeric integration strategy on finite domains
92#ifdef R__HAS_MATHMORE
93 RooAbsReal::defaultIntegratorConfig()->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D");
94
95 // 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
96 // ---------------------------------------------------------------------------------------
97
98 // Adjust maximum number of steps of RooIntegrator1D in the global default configuration
99 RooAbsReal::defaultIntegratorConfig()->getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 30);
100
101 // Example of how to change the parameters of a numeric integrator
102 // (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
103 // and RooCategories holding parameters with a finite set of options)
104 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg", 50);
105 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method", "15Points");
106
107 // Example of how to print set of possible values for "method" category
108 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method")->Print("v");
109#endif
110
111}
#define DEBUG
#define e(i)
Definition RSha256.hxx:103
double Double_t
Definition RtypesCore.h:59
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:231
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...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
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:91
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 ...
virtual Bool_t setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
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.
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:39
RooCmdArg Topic(Int_t topic)
RooCmdArg NumIntConfig(const RooNumIntConfig &cfg)
RooConstVar & RooConst(Double_t val)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...