Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf510_wsnamedsets.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #510
6##
7## Working with named parameter sets and parameter snapshots in
8## workspaces
9##
10## \macro_image
11## \macro_code
12## \macro_output
13##
14## \date February 2018
15## \authors Clemens Lange, Wouter Verkerke (C version)
16
17
18import ROOT
19
20
21def fillWorkspace(w):
22 # Create model
23 # -----------------------
24
25 # Declare observable x
26 x = ROOT.RooRealVar("x", "x", 0, 10)
27
28 # Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
29 # their parameters
30 mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
31 sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
32 sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
33
34 sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
35 sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
36
37 # Build Chebychev polynomial p.d.f.
38 a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
39 a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
40 bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
41
42 # Sum the signal components into a composite signal p.d.f.
43 sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
44 sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
45
46 # Sum the composite signal and background
47 bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
48 model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
49
50 # Import model into p.d.f.
51 w.Import(model)
52
53 # Encode definition of parameters in workspace
54 # ---------------------------------------------------------------------------------------
55
56 # Define named sets "parameters" and "observables", list which variables should be considered
57 # parameters and observables by the users convention
58 #
59 # Variables appearing in sets _must_ live in the workspace already, the autoImport flag
60 # of defineSet must be set to import them on the fly. Named sets contain only references
61 # to the original variables, the value of observables in named sets already
62 # reflect their 'current' value
63 params = model.getParameters({x})
64 w.defineSet("parameters", params)
65 w.defineSet("observables", {x})
66
67 # Encode reference value for parameters in workspace
68 # ---------------------------------------------------------------------------------------------------
69
70 # Define a parameter 'snapshot' in the p.d.f.
71 # Unlike a named set, parameter snapshot stores an independent set of values for
72 # a given set of variables in the workspace. The values can be stored and reloaded
73 # into the workspace variable objects using the loadSnapshot() and saveSnapshot()
74 # methods. A snapshot saves the value of each variable, errors that are stored
75 # with it as well as the 'Constant' flag that is used in fits to determine if a
76 # parameter is kept fixed or not.
77
78 # Do a dummy fit to a (supposedly) reference dataset here and store the results
79 # of that fit into a snapshot
80 refData = model.generate({x}, 10000)
81 model.fitTo(refData, PrintLevel=-1)
82
83 # The kTRUE flag imports the values of the objects in (*params) into the workspace
84 # If not set, present values of the workspace parameters objects are stored
85 w.saveSnapshot("reference_fit", params, True)
86
87 # Make another fit with the signal componentforced to zero
88 # and save those parameters too
89
90 bkgfrac.setVal(1)
91 bkgfrac.setConstant(True)
92 bkgfrac.removeError()
93 model.fitTo(refData, PrintLevel=-1)
94
95 w.saveSnapshot("reference_fit_bkgonly", params, True)
96
97
98# Create model and dataset
99# -----------------------------------------------
100
101w = ROOT.RooWorkspace("w")
102fillWorkspace(w)
103
104# Exploit convention encoded in named set "parameters" and "observables"
105# to use workspace contents w/o need for introspected
106model = w["model"]
107
108# Generate data from p.d.f. in given observables
109data = model.generate(w.set("observables"), 1000)
110
111# Fit model to data
112model.fitTo(data, PrintLevel=-1)
113
114# Plot fitted model and data on frame of first (only) observable
115frame = (w.set("observables").first()).frame()
116data.plotOn(frame)
117model.plotOn(frame)
118
119# Overlay plot with model with reference parameters as stored in snapshots
120w.loadSnapshot("reference_fit")
121model.plotOn(frame, LineColor="r")
122w.loadSnapshot("reference_fit_bkgonly")
123model.plotOn(frame, LineColor="r", LineStyle="--")
124
125# Draw the frame on the canvas
126c = ROOT.TCanvas("rf510_wsnamedsets", "rf503_wsnamedsets", 600, 600)
127ROOT.gPad.SetLeftMargin(0.15)
128frame.GetYaxis().SetTitleOffset(1.4)
129frame.Draw()
130
131c.SaveAs("rf510_wsnamedsets.png")
132
133# Print workspace contents
134w.Print()
135
136# Workspace will remain in memory after macro finishes
137ROOT.gDirectory.Add(w)
Definition first.py:1