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