Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf901_numintconfig.py
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_code
7## \macro_output
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 pdf integrals the relative change criterium is the most important,
28# however for certain non-pdf 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
40x = ROOT.RooRealVar("x", "x", -10, 10)
41landau = ROOT.RooLandau("landau", "landau", x, 0.0, 0.1)
42
43# Disable analytic integration from demonstration purposes
44landau.forceNumInt(True)
45
46# Activate debug-level messages for topic integration to be able to follow
47# actions below
48ROOT.RooMsgService.instance().addStream(ROOT.RooFit.DEBUG, Topic=ROOT.RooFit.Integration)
49
50# Calculate integral over landau with default choice of numeric integrator
51intLandau = landau.createIntegral({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(ROOT.RooAbsReal.defaultIntegratorConfig())
61integratorGKNotExisting = customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")
62if integratorGKNotExisting:
63 print("WARNING: RooAdaptiveGaussKronrodIntegrator is not existing because ROOT is built without Mathmore support")
64
65# Calculate integral over landau with custom integral specification
66intLandau2 = landau.createIntegral({x}, NumIntConfig=customConfig)
67val2 = intLandau2.getVal()
68print(" [2] int_dx landau(x) = ", val2)
69
70# Adjusting default config for a specific pdf
71# -------------------------------------------------------------------------------------
72
73# Another possibility: associate custom numeric integration configuration
74# as default for object 'landau'
75landau.setIntegratorConfig(customConfig)
76
77# Calculate integral over landau custom numeric integrator specified as
78# object default
79intLandau3 = landau.createIntegral({x})
80val3 = intLandau3.getVal()
81print(" [3] int_dx landau(x) = ", val3)
82
83# Another possibility: Change global default for 1D numeric integration
84# strategy on finite domains
85if not integratorGKNotExisting:
86 ROOT.RooAbsReal.defaultIntegratorConfig().method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")
87
88 # Adjusting parameters of a specific technique
89 # ---------------------------------------------------------------------------------------
90
91 # Adjust maximum number of steps of ROOT.RooIntegrator1D in the global
92 # default configuration
93 ROOT.RooAbsReal.defaultIntegratorConfig().getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 30)
94
95 # Example of how to change the parameters of a numeric integrator
96 # (Each config section is a ROOT.RooArgSet with ROOT.RooRealVars holding real-valued parameters
97 # and ROOT.RooCategories holding parameters with a finite set of options)
98 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg", 50)
99 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method", "15Points")
100
101 # Example of how to print set of possible values for "method" category
102 customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method").Print("v")