Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf901_numintconfig.py File Reference

Namespaces

namespace  rf901_numintconfig
 

Detailed Description

View in nbviewer Open in SWAN
Numeric algorithm tuning: configuration and customization of how numeric (partial) integrals are executed

from __future__ import print_function
import ROOT
# Adjust global 1D integration precision
# ----------------------------------------------------------------------------
# Print current global default configuration for numeric integration
# strategies
ROOT.RooAbsReal.defaultIntegratorConfig().Print("v")
# Example: Change global precision for 1D integrals from 1e-7 to 1e-6
#
# The relative epsilon (change as fraction of current best integral estimate) and
# absolute epsilon (absolute change w.r.t last best integral estimate) can be specified
# separately. For most pdf integrals the relative change criterium is the most important,
# however for certain non-pdf functions that integrate out to zero a separate absolute
# change criterium is necessary to declare convergence of the integral
#
# NB: ROOT.This change is for illustration only. In general the precision should be at least 1e-7
# for normalization integrals for MINUIT to succeed.
#
ROOT.RooAbsReal.defaultIntegratorConfig().setEpsAbs(1e-6)
ROOT.RooAbsReal.defaultIntegratorConfig().setEpsRel(1e-6)
# 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
# ------------------------------------------------------------------
x = ROOT.RooRealVar("x", "x", -10, 10)
landau = ROOT.RooLandau("landau", "landau", x, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(0.1))
# Disable analytic integration from demonstration purposes
landau.forceNumInt(True)
# Activate debug-level messages for topic integration to be able to follow
# actions below
ROOT.RooMsgService.instance().addStream(ROOT.RooFit.DEBUG, Topic=ROOT.RooFit.Integration)
# Calculate integral over landau with default choice of numeric integrator
intLandau = landau.createIntegral({x})
val = intLandau.getVal()
print(" [1] int_dx landau(x) = ", val) # setprecision(15)
# Same with custom configuration
# -----------------------------------------------------------
# Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
# for closed 1D integrals
customConfig = ROOT.RooNumIntConfig(ROOT.RooAbsReal.defaultIntegratorConfig())
integratorGKNotExisting = customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")
if integratorGKNotExisting:
print("WARNING: RooAdaptiveGaussKronrodIntegrator is not existing because ROOT is built without Mathmore support")
# Calculate integral over landau with custom integral specification
intLandau2 = landau.createIntegral({x}, NumIntConfig=customConfig)
val2 = intLandau2.getVal()
print(" [2] int_dx landau(x) = ", val2)
# Adjusting default config for a specific pdf
# -------------------------------------------------------------------------------------
# Another possibility: associate custom numeric integration configuration
# as default for object 'landau'
landau.setIntegratorConfig(customConfig)
# Calculate integral over landau custom numeric integrator specified as
# object default
intLandau3 = landau.createIntegral({x})
val3 = intLandau3.getVal()
print(" [3] int_dx landau(x) = ", val3)
# Another possibility: Change global default for 1D numeric integration
# strategy on finite domains
if not integratorGKNotExisting:
ROOT.RooAbsReal.defaultIntegratorConfig().method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")
# Adjusting parameters of a speficic technique
# ---------------------------------------------------------------------------------------
# Adjust maximum number of steps of ROOT.RooIntegrator1D in the global
# default configuration
ROOT.RooAbsReal.defaultIntegratorConfig().getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 30)
# Example of how to change the parameters of a numeric integrator
# (Each config section is a ROOT.RooArgSet with ROOT.RooRealVars holding real-valued parameters
# and ROOT.RooCategories holding parameters with a finite set of options)
customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setRealValue("maxSeg", 50)
customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").setCatLabel("method", "15Points")
# Example of how to print set of possible values for "method" category
customConfig.getConfigSection("RooAdaptiveGaussKronrodIntegrator1D").find("method").Print("v")
Requested precision: 1e-07 absolute, 1e-07 relative
1-D integration method: RooIntegrator1D (RooImproperIntegrator1D if open-ended)
2-D integration method: RooAdaptiveIntegratorND (N/A if open-ended)
N-D integration method: RooAdaptiveIntegratorND (N/A if open-ended)
Available integration methods:
*** RooBinIntegrator ***
Capabilities: [1-D] [2-D] [N-D]
Configuration:
1) numBins = 100
*** RooIntegrator1D ***
Capabilities: [1-D]
Configuration:
1) sumRule = Trapezoid(idx = 0)
2) extrapolation = Wynn-Epsilon(idx = 1)
3) maxSteps = 20
4) minSteps = 999
5) fixSteps = 0
*** RooIntegrator2D ***
Capabilities: [2-D]
Configuration:
(Depends on 'RooIntegrator1D')
*** RooSegmentedIntegrator1D ***
Capabilities: [1-D]
Configuration:
1) numSeg = 3
(Depends on 'RooIntegrator1D')
*** RooSegmentedIntegrator2D ***
Capabilities: [2-D]
Configuration:
(Depends on 'RooSegmentedIntegrator1D')
*** RooImproperIntegrator1D ***
Capabilities: [1-D] [OpenEnded]
Configuration:
(Depends on 'RooIntegrator1D')
*** RooMCIntegrator ***
Capabilities: [1-D] [2-D] [N-D]
Configuration:
1) samplingMode = Importance(idx = 0)
2) genType = QuasiRandom(idx = 0)
3) verbose = false(idx = 0)
4) alpha = 1.5
5) nRefineIter = 5
6) nRefinePerDim = 1000
7) nIntPerDim = 5000
*** RooAdaptiveIntegratorND ***
Capabilities: [2-D] [N-D]
Configuration:
1) maxEval2D = 100000
2) maxEval3D = 1e+06
3) maxEvalND = 1e+07
4) maxWarn = 5
*** RooAdaptiveGaussKronrodIntegrator1D ***
Capabilities: [1-D] [OpenEnded]
Configuration:
1) maxSeg = 100
2) method = 21Points(idx = 2)
*** RooGaussKronrodIntegrator1D ***
Capabilities: [1-D] [OpenEnded]
Configuration:
[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent
[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#3] INFO:Integration -- landau: Observables (x) are numerically integrated
[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent
[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#3] INFO:Integration -- landau: Observables (x) are numerically integrated
[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooAdaptiveGaussKronrodIntegrator1D to calculate Int(x)
[#3] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#3] DEBUG:Integration -- landau: Adding observable x as shape dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0 as value dependent
[#3] DEBUG:Integration -- landau: Adding parameter 0.1 as value dependent
[#3] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#3] INFO:Integration -- landau: Observables (x) are numerically integrated
[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooAdaptiveGaussKronrodIntegrator1D to calculate Int(x)
--- RooAbsArg ---
Value State: clean
Shape State: clean
Attributes:
Address: 0x133fc50
Clients:
Servers:
Proxies:
--- RooAbsCategory ---
Value = 1 "15Points)
Possible states:
[1] int_dx landau(x) = 0.09896533620544187
[2] int_dx landau(x) = 0.09895710292189497
[3] int_dx landau(x) = 0.09895710292189497
15Points 1
21Points 2
31Points 3
41Points 4
51Points 5
61Points 6
WynnEpsilon 0
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf901_numintconfig.py.