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