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