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