Logo ROOT   6.16/01
Reference Guide
rf201_composite.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #201
5## Composite p.d.f with signal and background component
6## pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
7##
8## \macro_code
9##
10## \date February 2018
11## \author Clemens Lange
12## \author Wouter Verkerke (C version)
13
14import ROOT
15
16# Setup component pdfs
17# ---------------------------------------
18
19# Declare observable x
20x = 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
24mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
25sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
26sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
27
28sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
29sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
30
31# Build Chebychev polynomial p.d.f.
32a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
33a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
34bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
35
36
37# Method 1 - Two RooAddPdfs
38# ------------------------------------------
39# Add signal components
40
41# Sum the signal components into a composite signal p.d.f.
42sig1frac = ROOT.RooRealVar(
43 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
44sig = ROOT.RooAddPdf("sig", "Signal", ROOT.RooArgList(
45 sig1, sig2), ROOT.RooArgList(sig1frac))
46
47# Add signal and background
48# ------------------------------------------------
49
50# Sum the composite signal and background
51bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
52model = ROOT.RooAddPdf(
53 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(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(ROOT.RooArgSet(x), 1000)
60
61# Fit model to data
62model.fitTo(data)
63
64# Plot data and PDF overlaid
65xframe = x.frame(ROOT.RooFit.Title(
66 "Example of composite pdf=(sig1+sig2)+bkg"))
67data.plotOn(xframe)
68model.plotOn(xframe)
69
70# Overlay the background component of model with a dashed line
71ras_bkg = ROOT.RooArgSet(bkg)
72model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg),
73 ROOT.RooFit.LineStyle(ROOT.kDashed))
74
75# Overlay the background+sig2 components of model with a dotted line
76ras_bkg_sig2 = ROOT.RooArgSet(bkg, sig2)
77model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg_sig2),
78 ROOT.RooFit.LineStyle(ROOT.kDotted))
79
80# Print structure of composite p.d.f.
81model.Print("t")
82
83# Method 2 - One RooAddPdf with recursive fractions
84# ---------------------------------------------------
85
86# Construct sum of models on one go using recursive fraction interpretations
87#
88# model2 = bkg + (sig1 + sig2)
89#
90model2 = ROOT.RooAddPdf("model", "g1+g2+a", ROOT.RooArgList(bkg,
91 sig1, sig2), ROOT.RooArgList(bkgfrac, sig1frac), ROOT.kTRUE)
92
93# NB: Each coefficient is interpreted as the fraction of the
94# left-hand component of the i-th recursive sum, i.e.
95#
96# sum4 = A + ( B + ( C + D) with fraction fA, and fC expands to
97#
98# sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
99
100# Plot recursive addition model
101# ---------------------------------------------------------
102model2.plotOn(xframe, ROOT.RooFit.LineColor(ROOT.kRed),
103 ROOT.RooFit.LineStyle(ROOT.kDashed))
104model2.plotOn(xframe, ROOT.RooFit.Components(ras_bkg_sig2),
105 ROOT.RooFit.LineColor(ROOT.kRed), ROOT.RooFit.LineStyle(ROOT.kDashed))
106model2.Print("t")
107
108# Draw the frame on the canvas
109c = ROOT.TCanvas("rf201_composite", "rf201_composite", 600, 600)
110ROOT.gPad.SetLeftMargin(0.15)
111xframe.GetYaxis().SetTitleOffset(1.4)
112xframe.Draw()
113
114c.SaveAs("rf201_composite.png")