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