Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf204a_extendedLikelihood.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Extended maximum likelihood fit in multiple ranges.
5##
6## \macro_code
7##
8## \date March 2021
9## \authors Harshal Shende, Stephan Hageboeck (C++ version)
10
11import ROOT
12
13
14# Setup component pdfs
15# ---------------------
16
17# Declare observable x
18x = ROOT.RooRealVar("x", "x", 0, 11)
19
20# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and 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
38# Extend the pdfs
39# -----------------------------
40
41# Define signal range in which events counts are to be defined
42x.setRange("signalRange", 4, 6)
43
44# Associated nsig/nbkg as expected number of events with sig/bkg _in_the_range_ "signalRange"
45nsig = ROOT.RooRealVar("nsig", "number of signal events in signalRange", 500, 0.0, 10000)
46nbkg = ROOT.RooRealVar("nbkg", "number of background events in signalRange", 500, 0, 10000)
47
48# Use AddPdf to extend the model. Giving as many coefficients as pdfs switches on extension.
49model = ROOT.RooAddPdf("model", "(g1+g2)+a", [bkg, sig], [nbkg, nsig])
50
51# Sample data, fit model
52# -------------------------------------------
53
54# Generate 1000 events from model so that nsig,nbkg come out to numbers <<500 in fit
55data = model.generate(x, 1000)
56
57canv = ROOT.TCanvas("Canvas", "Canvas", 1500, 600)
58canv.Divide(3, 1)
59
60# Fit full range
61# -------------------------------------------
62
63# Perform unbinned ML fit to data, full range
64
65# IMPORTANT:
66# The model needs to be copied when fitting with different ranges because
67# the interpretation of the coefficients is tied to the fit range
68# that's used in the first fit
69canv.cd(1)
70
71model1 = ROOT.RooAddPdf(model)
72r = model1.fitTo(data, Save=True)
73r.Print()
74
75frame = x.frame(Title="Full range fitted")
76data.plotOn(frame)
77model1.plotOn(frame, VisualizeError=r)
78model1.plotOn(frame)
79model1.paramOn(frame)
80frame.Draw()
81
82
83# Fit in two regions
84# -------------------------------------------
85
86canv.cd(2)
87x.setRange("left", 0.0, 4.0)
88x.setRange("right", 6.0, 10.0)
89
90model2 = ROOT.RooAddPdf(model)
91r2 = model2.fitTo(data, Range="left,right", Save=True)
92r2.Print()
93
94frame2 = x.frame(Title="Fit in left/right sideband")
95data.plotOn(frame2)
96model2.plotOn(frame2, VisualizeError=r2)
97model2.plotOn(frame2)
98model2.paramOn(frame2)
99frame2.Draw()
100
101
102# Fit in one region
103# -------------------------------------------
104# Note how restricting the region to only the left tail increases
105# the fit uncertainty
106
107canv.cd(3)
108x.setRange("leftToMiddle", 0.0, 5.0)
109
110model3 = ROOT.RooAddPdf(model)
111r3 = model3.fitTo(data, Range="leftToMiddle", Save=True)
112r3.Print()
113
114frame3 = x.frame(Title="Fit from left to middle")
115data.plotOn(frame3)
116model3.plotOn(frame3, VisualizeError=r3)
117model3.plotOn(frame3)
118model3.paramOn(frame3)
119frame3.Draw()
120
121canv.Draw()
122
123canv.SaveAs("rf204a_extendedLikelihood.png")