Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf205_compplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: options for plotting components of composite pdfs.
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Set up composite pdf
14# --------------------------------------
15
16# Declare observable x
17x = ROOT.RooRealVar("x", "x", 0, 10)
18
19# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
20# their parameters
21mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
22sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
23sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
24sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
25sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
26
27# Sum the signal components into a composite signal pdf
28sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
29sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
30
31# Build Chebychev polynomial pdf
32a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
33a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
34bkg1 = ROOT.RooChebychev("bkg1", "Background 1", x, [a0, a1])
35
36# Build expontential pdf
37alpha = ROOT.RooRealVar("alpha", "alpha", -1)
38bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
39
40# Sum the background components into a composite background pdf
41bkg1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
42bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [sig1frac])
43
44# Sum the composite signal and background
45bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
46model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
47
48# Set up basic plot with data and full pdf
49# ------------------------------------------------------------------------------
50
51# Generate a data sample of 1000 events in x from model
52data = model.generate({x}, 1000)
53
54# Plot data and complete PDF overlaid
55xframe = x.frame(Title="Component plotting of pdf=(sig1+sig2)+(bkg1+bkg2)")
56data.plotOn(xframe)
57model.plotOn(xframe)
58
59# Clone xframe for use below
60xframe2 = xframe.Clone("xframe2")
61
62# Make component by object reference
63# --------------------------------------------------------------------
64
65# Plot single background component specified by object reference
66ras_bkg = {bkg}
67model.plotOn(xframe, Components=ras_bkg, LineColor="r")
68
69# Plot single background component specified by object reference
70ras_bkg2 = {bkg2}
71model.plotOn(xframe, Components=ras_bkg2, LineStyle="--", LineColor="r")
72
73# Plot multiple background components specified by object reference
74# Note that specified components may occur at any level in object tree
75# (e.g bkg is component of 'model' and 'sig2' is component 'sig')
76ras_bkg_sig2 = {bkg, sig2}
77model.plotOn(xframe, Components=ras_bkg_sig2, LineStyle=":")
78
79# Make component by name/regexp
80# ------------------------------------------------------------
81
82# Plot single background component specified by name
83model.plotOn(xframe2, Components="bkg", LineColor="c")
84
85# Plot multiple background components specified by name
86model.plotOn(xframe2, Components="bkg1,sig2", LineStyle=":", LineColor="c")
87
88# Plot multiple background components specified by regular expression on
89# name
90model.plotOn(xframe2, Components="sig*", LineStyle="--", LineColor="c")
91
92# Plot multiple background components specified by multiple regular
93# expressions on name
94model.plotOn(xframe2, Invisible=True, Components="bkg1,sig*", LineStyle="--", LineColor="y")
95
96# Draw the frame on the canvas
97c = ROOT.TCanvas("rf205_compplot", "rf205_compplot", 800, 400)
98c.Divide(2)
99c.cd(1)
100ROOT.gPad.SetLeftMargin(0.15)
101xframe.GetYaxis().SetTitleOffset(1.4)
102xframe.Draw()
103c.cd(2)
104ROOT.gPad.SetLeftMargin(0.15)
105xframe2.GetYaxis().SetTitleOffset(1.4)
106xframe2.Draw()
107
108c.SaveAs("rf205_compplot.png")