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