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., 1.)
30a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
31bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
32
33# Sum the signal components into a composite signal pdf
34sig1frac = ROOT.RooRealVar(
35 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
36sig = ROOT.RooAddPdf(
37 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
38
39# Method 1 - Construct extended composite model
40# -------------------------------------------------------------------
41
42# Sum the composite signal and background into an extended pdf
43# nsig*sig+nbkg*bkg
44nsig = ROOT.RooRealVar("nsig", "number of signal events", 500, 0., 10000)
45nbkg = ROOT.RooRealVar(
46 "nbkg", "number of background events", 500, 0, 10000)
47model = ROOT.RooAddPdf(
48 "model",
49 "(g1+g2)+a",
50 ROOT.RooArgList(
51 bkg,
52 sig),
53 ROOT.RooArgList(
54 nbkg,
55 nsig))
56
57# Sample, fit and plot extended model
58# ---------------------------------------------------------------------
59
60# Generate a data sample of expected number events in x from model
61# = model.expectedEvents() = nsig+nbkg
62data = model.generate(ROOT.RooArgSet(x))
63
64# Fit model to data, ML term automatically included
65model.fitTo(data)
66
67# Plot data and PDF overlaid, expected number of events for pdf projection normalization
68# rather than observed number of events (==data.numEntries())
69xframe = x.frame(ROOT.RooFit.Title("extended ML fit example"))
70data.plotOn(xframe)
71model.plotOn(xframe, ROOT.RooFit.Normalization(
72 1.0, ROOT.RooAbsReal.RelativeExpected))
73
74# Overlay the background component of model with a dashed line
75ras_bkg = ROOT.RooArgSet(bkg)
76model.plotOn(
77 xframe, ROOT.RooFit.Components(ras_bkg), ROOT.RooFit.LineStyle(
78 ROOT.kDashed), ROOT.RooFit.Normalization(
79 1.0, ROOT.RooAbsReal.RelativeExpected))
80
81# Overlay the background+sig2 components of model with a dotted line
82ras_bkg_sig2 = ROOT.RooArgSet(bkg, sig2)
83model.plotOn(
84 xframe, ROOT.RooFit.Components(ras_bkg_sig2), ROOT.RooFit.LineStyle(
85 ROOT.kDotted), ROOT.RooFit.Normalization(
86 1.0, ROOT.RooAbsReal.RelativeExpected))
87
88# Print structure of composite pdf
89model.Print("t")
90
91
92# Method 2 - Construct extended components first
93# ---------------------------------------------------------------------
94
95# Associated nsig/nbkg as expected number of events with sig/bkg
96esig = ROOT.RooExtendPdf("esig", "extended signal pdf", sig, nsig)
97ebkg = ROOT.RooExtendPdf("ebkg", "extended background pdf", bkg, nbkg)
98
99# Sum extended components without coefs
100# -------------------------------------------------------------------------
101
102# Construct sum of two extended pdf (no coefficients required)
103model2 = ROOT.RooAddPdf("model2", "(g1+g2)+a", ROOT.RooArgList(ebkg, esig))
104
105# Draw the frame on the canvas
106c = ROOT.TCanvas("rf202_extendedmlfit", "rf202_extendedmlfit", 600, 600)
107ROOT.gPad.SetLeftMargin(0.15)
108xframe.GetYaxis().SetTitleOffset(1.4)
109xframe.Draw()
110
111c.SaveAs("rf202_extendedmlfit.png")