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