Logo ROOT   6.12/07
Reference Guide
rf901_numintconfig.C File Reference

Detailed Description

View in nbviewer Open in SWAN 'NUMERIC ALGORITHM TUNING' RooFit tutorial macro #901

Configuration and customization of how numeric (partial) integrals are executed

0.0293860435486
4.80892395973
Processing /mnt/build/workspace/root-makedoc-v612/rootspi/rdoc/src/v6-12-00-patches/tutorials/roofit/rf901_numintconfig.C...
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
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
*** RooAdaptiveGaussKronrodIntegrator1D ***
Capabilities: [1-D] [OpenEnded]
Configuration:
1) maxSeg = 100
2) method = 21Points(idx = 2)
*** RooGaussKronrodIntegrator1D ***
Capabilities: [1-D] [OpenEnded]
Configuration:
*** RooAdaptiveIntegratorND ***
Capabilities: [2-D] [N-D]
Configuration:
1) maxEval2D = 100000
2) maxEval3D = 1e+06
3) maxEvalND = 1e+07
4) maxWarn = 5
[#2] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#2] DEBUG:Integration -- landau: Adding observable x of server x as shape dependent
[#2] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#2] INFO:Integration -- landau: Observables (x) are numerically integrated
[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[1] int_dx landau(x) = 0.0989653362054419
[#2] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#2] DEBUG:Integration -- landau: Adding observable x of server x as shape dependent
[#2] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#2] INFO:Integration -- landau: Observables (x) are numerically integrated
[#1] INFO:NumericIntegration -- RooRealIntegral::init(landau_Int[x]) using numeric integrator RooAdaptiveGaussKronrodIntegrator1D to calculate Int(x)
[2] int_dx landau(x) = 0.098957102921895
[#2] INFO:Integration -- RooRealIntegral::ctor(landau_Int[x]) Constructing integral of function landau over observables(x) with normalization () with range identifier <none>
[#2] DEBUG:Integration -- landau: Adding observable x of server x as shape dependent
[#2] INFO:Integration -- landau: Observable x is suitable for analytical integration (if supported by p.d.f)
[#2] 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] int_dx landau(x) = 0.098957102921895
--- RooAbsArg ---
Value State: clean
Shape State: clean
Attributes:
Address: 0x397bad0
Clients:
Servers:
Proxies:
--- RooAbsCategory ---
Value is "15Points" (1)
Has the following possible values:
WynnEpsilon = 0
15Points = 1
21Points = 2
31Points = 3
41Points = 4
51Points = 5
61Points = 6
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "RooLandau.h"
#include "RooArgSet.h"
#include <iomanip>
using namespace RooFit ;
void rf901_numintconfig()
{
// A d j u s t g l o b a l 1 D i n t e g r a t i o n p r e c i s i o n
// ----------------------------------------------------------------------------
// Print current global default configuration for numeric integration strategies
// 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 p.d.f integrals the relative change criterium is the most important,
// however for certain non-p.d.f functions that integrate out to zero a separate absolute
// change criterium is necessary to declare convergence of the integral
//
// NB: This change is for illustration only. In general the precision should be at least 1e-7
// for normalization integrals for MINUIT to succeed.
//
// 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
// ------------------------------------------------------------------
// Construct p.d.f without support for analytical integrator for demonstration purposes
RooRealVar x("x","x",-10,10) ;
RooLandau landau("landau","landau",x,RooConst(0),RooConst(0.1)) ;
// Activate debug-level messages for topic integration to be able to follow actions below
// Calculate integral over landau with default choice of numeric integrator
RooAbsReal* intLandau = landau.createIntegral(x) ;
Double_t val = intLandau->getVal() ;
cout << " [1] int_dx landau(x) = " << setprecision(15) << val << endl ;
// S a m e w i t h c u s t o m c o n f i g u r a t i o n
// -----------------------------------------------------------
// Construct a custom configuration which uses the adaptive Gauss-Kronrod technique
// for closed 1D integrals
customConfig.method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
// Calculate integral over landau with custom integral specification
RooAbsReal* intLandau2 = landau.createIntegral(x,NumIntConfig(customConfig)) ;
Double_t val2 = intLandau2->getVal() ;
cout << " [2] int_dx landau(x) = " << val2 << endl ;
// A d j u s t i n g d e f a u l t c o n f i g f o r a s p e c i f i c p d f
// -------------------------------------------------------------------------------------
// 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
RooAbsReal* intLandau3 = landau.createIntegral(x) ;
Double_t val3 = intLandau3->getVal() ;
cout << " [3] int_dx landau(x) = " << val3 << endl ;
// Another possibility: Change global default for 1D numeric integration strategy on finite domains
RooAbsReal::defaultIntegratorConfig()->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
// A d j u s t i n g p a r a m e t e r s o f a s p e c i f i c t e c h n i q u e
// ---------------------------------------------------------------------------------------
// Adjust maximum number of steps of RooIntegrator1D in the global default configuration
RooAbsReal::defaultIntegratorConfig()->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",30) ;
// Example of how to change the parameters of a numeric integrator
// (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
// and 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") ;
}
Author
07/2008 - Wouter Verkerke

Definition in file rf901_numintconfig.C.