Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf509_wsinteractive.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Organization and simultaneous fits: easy interactive access to workspace contents - CINT
5## to CLING code migration
6##
7## \macro_image
8## \macro_code
9## \macro_output
10##
11## \date February 2018
12## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14import ROOT
15
16
17def fillWorkspace(w):
18 # Create pdf and fill workspace
19 # --------------------------------------------------------
20
21 # Declare observable x
22 x = ROOT.RooRealVar("x", "x", 0, 10)
23
24 # Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
25 # their parameters
26 mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
27 sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
28 sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
29
30 sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
31 sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
32
33 # Build Chebychev polynomial pdf
34 a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
35 a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
36 bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
37
38 # Sum the signal components into a composite signal pdf
39 sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
40 sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
41
42 # Sum the composite signal and background
43 bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
44 model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
45
46 w.Import(model)
47
48
49# Create and fill workspace
50# ------------------------------------------------
51
52
53# Create a workspace named 'w'
54# With CINT w could exports its contents to
55# a same-name C++ namespace in CINT 'namespace w'.
56# but self does not work anymore in CLING.
57# so self tutorial is an example on how to
58# change the code
59w = ROOT.RooWorkspace("w", True)
60
61# Fill workspace with pdf and data in a separate function
62fillWorkspace(w)
63
64# Print workspace contents
65w.Print()
66
67# self does not work anymore with CLING
68# use normal workspace functionality
69
70# Use workspace contents
71# ----------------------------------------------
72
73# Old syntax to use the name space prefix operator to access the workspace contents
74#
75# d = w.model.generate(w.x,1000)
76# r = w.model.fitTo(*d)
77
78# use normal workspace methods
79model = w["model"]
80x = w["x"]
81
82d = model.generate({x}, 1000)
83r = model.fitTo(d, PrintLevel=-1)
84
85# old syntax to access the variable x
86# frame = w.x.frame()
87
88frame = x.frame()
89d.plotOn(frame)
90
91# OLD syntax to omit x.
92# NB: The 'w.' prefix can be omitted if namespace w is imported in local namespace
93# in the usual C++ way
94#
95# using namespace w
96# model.plotOn(frame)
97# model.plotOn(frame, Components=bkg, LineStyle="--")
98
99# correct syntax
100bkg = w["bkg"]
101model.plotOn(frame)
102model.plotOn(frame, Components=bkg, LineStyle="--")
103
104# Draw the frame on the canvas
105c = ROOT.TCanvas("rf509_wsinteractive", "rf509_wsinteractive", 600, 600)
106ROOT.gPad.SetLeftMargin(0.15)
107frame.GetYaxis().SetTitleOffset(1.4)
108frame.Draw()
109
110c.SaveAs("rf509_wsinteractive.png")