Logo ROOT   6.16/01
Reference Guide
rf903_numintcache.py
Go to the documentation of this file.
1## \ingroup tutorial_roofit
2## \notebook
3##
4## 'NUMERIC ALGORITHM TUNING' RooFit tutorial macro #903
5##
6## Caching of slow numeric integrals and parameterizations of slow
7## numeric integrals
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13
14
15import sys
16import ROOT
17
18
19def getWorkspace(mode):
20 # Create, save or load workspace with pdf
21 # -----------------------------------------------------------------------------------
22 #
23 # Mode = 0 : Create workspace for plain running (no integral caching)
24 # Mode = 1 : Generate workspace with precalculated integral and store it on file
25 # Mode = 2 : Load previously stored workspace from file
26
27 w = ROOT.RooWorkspace()
28
29 if mode != 2:
30 # Create empty workspace workspace
31 w = ROOT.RooWorkspace("w", 1)
32
33 # Make a difficult to normalize p.d.f. in 3 dimensions that is
34 # integrated numerically.
35 w.factory(
36 "EXPR::model('1/((x-a)*(x-a)+0.01)+1/((y-a)*(y-a)+0.01)+1/((z-a)*(z-a)+0.01)',x[-1,1],y[-1,1],z[-1,1],a[-5,5])")
37
38 if mode == 1:
39 # Instruct model to precalculate normalization integral that integrate at least
40 # two dimensions numerically. In self specific case the integral value for
41 # all values of parameter 'a' are stored in a histogram and available for use
42 # in subsequent fitting and plotting operations (interpolation is
43 # applied)
44
45 # w.pdf("model").setNormValueCaching(3)
46 w.pdf("model").setStringAttribute("CACHEPARMINT", "x:y:z")
47
48 # Evaluate p.d.f. once to trigger filling of cache
49 normSet = ROOT.RooArgSet(w.var("x"), w.var("y"), w.var("z"))
50 w.pdf("model").getVal(normSet)
51 w.writeToFile("rf903_numintcache.root")
52
53 if (mode == 2):
54 # Load preexisting workspace from file in mode==2
55 f = ROOT.TFile("rf903_numintcache.root")
56 w = f.Get("w")
57
58 # Return created or loaded workspace
59 return w
60
61
62mode=0
63# Mode = 0 : Run plain fit (slow)
64# Mode = 1 : Generate workspace with precalculated integral and store it on file (prepare for accelerated running)
65# Mode = 2 : Run fit from previously stored workspace including cached
66# integrals (fast, run in mode=1 first)
67
68# Create, save or load workspace with pdf
69# -----------------------------------------------------------------------------------
70
71# Make/load workspace, here in mode 1
72w = getWorkspace(mode)
73if mode == 1:
74 # Show workspace that was created
75 w.Print()
76
77 # Show plot of cached integral values
78 hhcache = w.expensiveObjectCache().getObj(1)
79 if (hhcache):
80 ROOT.TCanvas("rf903_numintcache", "rf903_numintcache", 600, 600)
81 hhcache.createHistogram("a").Draw()
82 else:
83 ROOT.RooFit.Error("rf903_numintcache",
84 "Cached histogram is not existing in workspace")
85 sys.exit()
86
87# Use pdf from workspace for generation and fitting
88# -----------------------------------------------------------------------------------
89
90# ROOT.This is always slow (need to find maximum function value
91# empirically in 3D space)
92d = w.pdf("model").generate(ROOT.RooArgSet(w.var("x"), w.var("y"), w.var("z")), 1000)
93
94# ROOT.This is slow in mode 0, fast in mode 1
95w.pdf("model").fitTo(d, ROOT.RooFit.Verbose(ROOT.kTRUE), ROOT.RooFit.Timer(ROOT.kTRUE))
96
97# Projection on x (always slow as 2D integral over Y, at fitted value of a
98# is not cached)
99framex = w.var("x").frame(ROOT.RooFit.Title("Projection of 3D model on X"))
100d.plotOn(framex)
101w.pdf("model").plotOn(framex)
102
103# Draw x projection on canvas
104c = ROOT.TCanvas("rf903_numintcache", "rf903_numintcache", 600, 600)
105framex.Draw()
106
107c.SaveAs("rf903_numintcache.png")
108
109# Make workspace available on command line after macro finishes
110ROOT.gDirectory.Add(w)
th1 Draw()