Logo ROOT   6.16/01
Reference Guide
rf202_extendedmlfit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #202
5## Setting up an extended maximum likelihood fit
6##
7## \macro_code
8##
9## \date February 2018
10## \author Clemens Lange
11## \author 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 p.d.f.
31a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
32a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
33bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
34
35# Sum the signal components into a composite signal p.d.f.
36sig1frac = ROOT.RooRealVar(
37 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
38sig = ROOT.RooAddPdf(
39 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
40
41# Method 1 - Construct extended composite model
42# -------------------------------------------------------------------
43
44# Sum the composite signal and background into an extended pdf
45# nsig*sig+nbkg*bkg
46nsig = ROOT.RooRealVar("nsig", "number of signal events", 500, 0., 10000)
47nbkg = ROOT.RooRealVar(
48 "nbkg", "number of background events", 500, 0, 10000)
49model = ROOT.RooAddPdf(
50 "model", "(g1+g2)+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(nbkg, nsig))
51
52# Sample, fit and plot extended model
53# ---------------------------------------------------------------------
54
55# Generate a data sample of expected number events in x from model
56# = model.expectedEvents() = nsig+nbkg
57data = model.generate(ROOT.RooArgSet(x))
58
59# Fit model to data, ML term automatically included
60model.fitTo(data)
61
62# Plot data and PDF overlaid, expected number of events for p.d.f projection normalization
63# rather than observed number of events (==data.numEntries())
64xframe = x.frame(ROOT.RooFit.Title("extended ML fit example"))
65data.plotOn(xframe)
66model.plotOn(xframe, ROOT.RooFit.Normalization(
67 1.0, ROOT.RooAbsReal.RelativeExpected))
68
69# Overlay the background component of model with a dashed line
70ras_bkg = ROOT.RooArgSet(bkg)
71model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg), ROOT.RooFit.LineStyle(ROOT.kDashed),
72 ROOT.RooFit.Normalization(1.0, ROOT.RooAbsReal.RelativeExpected))
73
74# Overlay the background+sig2 components of model with a dotted line
75ras_bkg_sig2 = ROOT.RooArgSet(bkg, sig2)
76model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg_sig2), ROOT.RooFit.LineStyle(
77 ROOT.kDotted), ROOT.RooFit.Normalization(1.0, ROOT.RooAbsReal.RelativeExpected))
78
79# Print structure of composite p.d.f.
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 p.d.f", sig, nsig)
88ebkg = ROOT.RooExtendPdf("ebkg", "extended background p.d.f", bkg, nbkg)
89
90# Sum extended components without coefs
91# -------------------------------------------------------------------------
92
93# Construct sum of two extended p.d.f. (no coefficients required)
94model2 = ROOT.RooAddPdf("model2", "(g1+g2)+a", ROOT.RooArgList(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")