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