Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf902_numgenconfig.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 MC sampling algorithms
5## on specific pdfs are executed
6##
7## \macro_code
8## \macro_output
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Adjust global MC sampling strategy
17# ------------------------------------------------------------------
18
19# Example pdf for use below
20x = ROOT.RooRealVar("x", "x", 0, 10)
21model = ROOT.RooChebychev("model", "model", x, [0.0, 0.5, -0.1])
22
23# Change global strategy for 1D sampling problems without conditional observable
24# (1st kFALSE) and without discrete observable (2nd kFALSE) from ROOT.RooFoamGenerator,
25# ( an interface to the ROOT.TFoam MC generator with adaptive subdivisioning strategy ) to ROOT.RooAcceptReject,
26# a plain accept/reject sampling algorithm [ ROOT.RooFit default before
27# ROOT 5.23/04 ]
28ROOT.RooAbsPdf.defaultGeneratorConfig().method1D(False, False).setLabel("RooAcceptReject")
29
30# Generate 10Kevt using ROOT.RooAcceptReject
31data_ar = model.generate({x}, 10000, Verbose=True)
32data_ar.Print()
33
34# Adjusting default config for a specific pdf
35# -------------------------------------------------------------------------------------
36
37# Another possibility: associate custom MC sampling configuration as default for object 'model'
38# The kTRUE argument will install a clone of the default configuration as specialized configuration
39# for self model if none existed so far
40model.specialGeneratorConfig(True).method1D(False, False).setLabel("RooFoamGenerator")
41
42# Adjusting parameters of a specific technique
43# ---------------------------------------------------------------------------------------
44
45# Adjust maximum number of steps of ROOT.RooIntegrator1D in the global
46# default configuration
47ROOT.RooAbsPdf.defaultGeneratorConfig().getConfigSection("RooAcceptReject").setRealValue("nTrial1D", 2000)
48
49# Example of how to change the parameters of a numeric integrator
50# (Each config section is a ROOT.RooArgSet with ROOT.RooRealVars holding real-valued parameters
51# and ROOT.RooCategories holding parameters with a finite set of options)
52model.specialGeneratorConfig().getConfigSection("RooFoamGenerator").setRealValue("chatLevel", 1)
53
54# Generate 10Kevt using ROOT.RooFoamGenerator (FOAM verbosity increased
55# with above chatLevel adjustment for illustration purposes)
56data_foam = model.generate({x}, 10000, Verbose=True)
57data_foam.Print()