Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf202_extendedmlfit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: setting up an extended maximum likelihood fit
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Set up component pdfs
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)
24
25sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
26sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
27
28# Build Chebychev polynomial pdf
29a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
30a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
31bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
32
33# Sum the signal components into a composite signal pdf
34sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
35sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
36
37# Method 1 - Construct extended composite model
38# -------------------------------------------------------------------
39
40# Sum the composite signal and background into an extended pdf
41# nsig*sig+nbkg*bkg
42nsig = ROOT.RooRealVar("nsig", "number of signal events", 500, 0.0, 10000)
43nbkg = ROOT.RooRealVar("nbkg", "number of background events", 500, 0, 10000)
44model = ROOT.RooAddPdf("model", "(g1+g2)+a", [bkg, sig], [nbkg, nsig])
45
46# Sample, fit and plot extended model
47# ---------------------------------------------------------------------
48
49# Generate a data sample of expected number events in x from model
50# = model.expectedEvents() = nsig+nbkg
51data = model.generate({x})
52
53# Fit model to data, ML term automatically included
54model.fitTo(data)
55
56# Plot data and PDF overlaid, expected number of events for pdf projection normalization
57# rather than observed number of events (==data.numEntries())
58xframe = x.frame(Title="extended ML fit example")
59data.plotOn(xframe)
60model.plotOn(xframe, Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected))
61
62# Overlay the background component of model with a dashed line
63model.plotOn(
64 xframe,
65 Components={bkg},
66 LineStyle=":",
67 Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected),
68)
69
70# Overlay the background+sig2 components of model with a dotted line
71ras_bkg_sig2 = {bkg, sig2}
72model.plotOn(
73 xframe,
74 Components=ras_bkg_sig2,
75 LineStyle=":",
76 Normalization=dict(scaleFactor=1.0, scaleType=ROOT.RooAbsReal.RelativeExpected),
77)
78
79# Print structure of composite pdf
80model.Print("t")
81
82
83# Method 2 - Construct extended components first
84# ---------------------------------------------------------------------
85
86# Associated nsig/nbkg as expected number of events with sig/bkg
87esig = ROOT.RooExtendPdf("esig", "extended signal pdf", sig, nsig)
88ebkg = ROOT.RooExtendPdf("ebkg", "extended background pdf", bkg, nbkg)
89
90# Sum extended components without coefs
91# -------------------------------------------------------------------------
92
93# Construct sum of two extended pdf (no coefficients required)
94model2 = ROOT.RooAddPdf("model2", "(g1+g2)+a", [ebkg, esig])
95
96# Draw the frame on the canvas
97c = ROOT.TCanvas("rf202_extendedmlfit", "rf202_extendedmlfit", 600, 600)
98ROOT.gPad.SetLeftMargin(0.15)
99xframe.GetYaxis().SetTitleOffset(1.4)
100xframe.Draw()
101
102c.SaveAs("rf202_extendedmlfit.png")