Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf201_composite.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: composite pdf with signal and background component
5##
6## ```
7## pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
8## ```
9##
10## \macro_image
11## \macro_code
12## \macro_output
13##
14## \date February 2018
15## \authors Clemens Lange, Wouter Verkerke (C++ version)
16
17import ROOT
18
19# Setup component pdfs
20# ---------------------------------------
21
22# Declare observable x
23x = ROOT.RooRealVar("x", "x", 0, 10)
24
25# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
26# their parameters
27mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
28sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
29sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
30
31sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
32sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
33
34# Build Chebychev polynomial pdf
35a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
36a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
37bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
38
39
40# Method 1 - Two RooAddPdfs
41# ------------------------------------------
42# Add signal components
43
44# Sum the signal components into a composite signal pdf
45sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
46sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
47
48# Add signal and background
49# ------------------------------------------------
50
51# Sum the composite signal and background
52bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
53model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
54
55# Sample, fit and plot model
56# ---------------------------------------------------
57
58# Generate a data sample of 1000 events in x from model
59data = model.generate({x}, 1000)
60
61# Fit model to data
62model.fitTo(data, PrintLevel=-1)
63
64# Plot data and PDF overlaid
65xframe = x.frame(Title="Example of composite pdf=(sig1+sig2)+bkg")
66data.plotOn(xframe)
67model.plotOn(xframe)
68
69# Overlay the background component of model with a dashed line
70model.plotOn(xframe, Components={bkg}, LineStyle="--")
71
72# Overlay the background+sig2 components of model with a dotted line
73model.plotOn(xframe, Components={bkg, sig2}, LineStyle=":")
74
75# Print structure of composite pdf
76model.Print("t")
77
78# Method 2 - One RooAddPdf with recursive fractions
79# ---------------------------------------------------
80
81# Construct sum of models on one go using recursive fraction interpretations
82#
83# model2 = bkg + (sig1 + sig2)
84#
85model2 = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig1, sig2], [bkgfrac, sig1frac], True)
86
87# NB: Each coefficient is interpreted as the fraction of the
88# left-hand component of the i-th recursive sum, i.e.
89#
90# sum4 = A + ( B + ( C + D) with fraction fA, and fC expands to
91#
92# sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
93
94# Plot recursive addition model
95# ---------------------------------------------------------
96model2.plotOn(xframe, LineColor="r", LineStyle="--")
97model2.plotOn(xframe, Components={bkg, sig2}, LineColor="r", LineStyle="--")
98model2.Print("t")
99
100# Draw the frame on the canvas
101c = ROOT.TCanvas("rf201_composite", "rf201_composite", 600, 600)
102ROOT.gPad.SetLeftMargin(0.15)
103xframe.GetYaxis().SetTitleOffset(1.4)
104xframe.Draw()
105
106c.SaveAs("rf201_composite.png")