Logo ROOT   6.16/01
Reference Guide
rf610_visualerror.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #610
6##
7## Visualization of errors from a covariance matrix
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15
16import ROOT
17
18
19# Setup example fit
20# ---------------------------------------
21
22# Create sum of two Gaussians p.d.f. with factory
23x = ROOT.RooRealVar("x", "x", -10, 10)
24
25m = ROOT.RooRealVar("m", "m", 0, -10, 10)
26s = ROOT.RooRealVar("s", "s", 2, 1, 50)
27sig = ROOT.RooGaussian("sig", "sig", x, m, s)
28
29m2 = ROOT.RooRealVar("m2", "m2", -1, -10, 10)
30s2 = ROOT.RooRealVar("s2", "s2", 6, 1, 50)
31bkg = ROOT.RooGaussian("bkg", "bkg", x, m2, s2)
32
33fsig = ROOT.RooRealVar("fsig", "fsig", 0.33, 0, 1)
34model = ROOT.RooAddPdf("model", "model", ROOT.RooArgList(
35 sig, bkg), ROOT.RooArgList(fsig))
36
37# Create binned dataset
38x.setBins(25)
39d = model.generateBinned(ROOT.RooArgSet(x), 1000)
40
41# Perform fit and save fit result
42r = model.fitTo(d, ROOT.RooFit.Save())
43
44# Visualize fit error
45# -------------------------------------
46
47# Make plot frame
48frame = x.frame(ROOT.RooFit.Bins(40), ROOT.RooFit.Title(
49 "P.d.f with visualized 1-sigma error band"))
50d.plotOn(frame)
51
52# Visualize 1-sigma error encoded in fit result 'r' as orange band using linear error propagation
53# ROOT.This results in an error band that is by construction symmetric
54#
55# The linear error is calculated as
56# error(x) = Z* F_a(x) * Corr(a,a') F_a'(x)
57#
58# where F_a(x) = [ f(x,a+da) - f(x,a-da) ] / 2,
59#
60# with f(x) = the plotted curve
61# 'da' = error taken from the fit result
62# Corr(a,a') = the correlation matrix from the fit result
63# Z = requested significance 'Z sigma band'
64#
65# The linear method is fast (required 2*N evaluations of the curve, N is the number of parameters),
66# but may not be accurate in the presence of strong correlations (~>0.9) and at Z>2 due to linear and
67# Gaussian approximations made
68#
69model.plotOn(frame, ROOT.RooFit.VisualizeError(
70 r, 1), ROOT.RooFit.FillColor(ROOT.kOrange))
71
72# Calculate error using sampling method and visualize as dashed red line.
73#
74# In self method a number of curves is calculated with variations of the parameter values, sampled
75# from a multi-variate Gaussian p.d.f. that is constructed from the fit results covariance matrix.
76# The error(x) is determined by calculating a central interval that capture N% of the variations
77# for each valye of x, N% is controlled by Z (i.e. Z=1 gives N=68%). The number of sampling curves
78# is chosen to be such that at least 100 curves are expected to be outside the N% interval, is minimally
79# 100 (e.g. Z=1.Ncurve=356, Z=2.Ncurve=2156)) Intervals from the sampling method can be asymmetric,
80# and may perform better in the presence of strong correlations, may take
81# (much) longer to calculate
82model.plotOn(frame, ROOT.RooFit.VisualizeError(r, 1, ROOT.kFALSE), ROOT.RooFit.DrawOption(
83 "L"), ROOT.RooFit.LineWidth(2), ROOT.RooFit.LineColor(ROOT.kRed))
84
85# Perform the same type of error visualization on the background component only.
86# The VisualizeError() option can generally applied to _any_ kind of
87# plot (components, asymmetries, etc..)
88model.plotOn(frame, ROOT.RooFit.VisualizeError(r, 1),
89 ROOT.RooFit.FillColor(ROOT.kOrange), ROOT.RooFit.Components("bkg"))
90model.plotOn(frame, ROOT.RooFit.VisualizeError(r, 1, ROOT.kFALSE), ROOT.RooFit.DrawOption("L"), ROOT.RooFit.LineWidth(
91 2), ROOT.RooFit.LineColor(ROOT.kRed), ROOT.RooFit.Components("bkg"), ROOT.RooFit.LineStyle(ROOT.kDashed))
92
93# Overlay central value
94model.plotOn(frame)
95model.plotOn(frame, ROOT.RooFit.Components("bkg"),
96 ROOT.RooFit.LineStyle(ROOT.kDashed))
97d.plotOn(frame)
98frame.SetMinimum(0)
99
100# Visualize partial fit error
101# ------------------------------------------------------
102
103# Make plot frame
104frame2 = x.frame(ROOT.RooFit.Bins(40), ROOT.RooFit.Title(
105 "Visualization of 2-sigma partial error from (m,m2)"))
106
107# Visualize partial error. For partial error visualization the covariance matrix is first reduced as follows
108# ___ -1
109# Vred = V22 = V11 - V12 * V22 * V21
110#
111# Where V11,V12,V21, represent a block decomposition of the covariance matrix into observables that
112# are propagated (labeled by index '1') and that are not propagated (labeled by index '2'), V22bar
113# is the Shur complement of V22, as shown above
114#
115# (Note that Vred is _not_ a simple sub-matrix of V)
116
117# Propagate partial error due to shape parameters (m,m2) using linear and
118# sampling method
119model.plotOn(frame2, ROOT.RooFit.VisualizeError(
120 r, ROOT.RooArgSet(m, m2), 2), ROOT.RooFit.FillColor(ROOT.kCyan))
121model.plotOn(frame2, ROOT.RooFit.Components("bkg"), ROOT.RooFit.VisualizeError(
122 r, ROOT.RooArgSet(m, m2), 2), ROOT.RooFit.FillColor(ROOT.kCyan))
123
124model.plotOn(frame2)
125model.plotOn(frame2, ROOT.RooFit.Components("bkg"),
126 ROOT.RooFit.LineStyle(ROOT.kDashed))
127frame2.SetMinimum(0)
128
129# Make plot frame
130frame3 = x.frame(ROOT.RooFit.Bins(40), ROOT.RooFit.Title(
131 "Visualization of 2-sigma partial error from (s,s2)"))
132
133# Propagate partial error due to yield parameter using linear and sampling
134# method
135model.plotOn(frame3, ROOT.RooFit.VisualizeError(
136 r, ROOT.RooArgSet(s, s2), 2), ROOT.RooFit.FillColor(ROOT.kGreen))
137model.plotOn(frame3, ROOT.RooFit.Components("bkg"), ROOT.RooFit.VisualizeError(
138 r, ROOT.RooArgSet(s, s2), 2), ROOT.RooFit.FillColor(ROOT.kGreen))
139
140model.plotOn(frame3)
141model.plotOn(frame3, ROOT.RooFit.Components("bkg"),
142 ROOT.RooFit.LineStyle(ROOT.kDashed))
143frame3.SetMinimum(0)
144
145# Make plot frame
146frame4 = x.frame(ROOT.RooFit.Bins(40), ROOT.RooFit.Title(
147 "Visualization of 2-sigma partial error from fsig"))
148
149# Propagate partial error due to yield parameter using linear and sampling
150# method
151model.plotOn(frame4, ROOT.RooFit.VisualizeError(
152 r, ROOT.RooArgSet(fsig), 2), ROOT.RooFit.FillColor(ROOT.kMagenta))
153model.plotOn(frame4, ROOT.RooFit.Components("bkg"), ROOT.RooFit.VisualizeError(
154 r, ROOT.RooArgSet(fsig), 2), ROOT.RooFit.FillColor(ROOT.kMagenta))
155
156model.plotOn(frame4)
157model.plotOn(frame4, ROOT.RooFit.Components("bkg"),
158 ROOT.RooFit.LineStyle(ROOT.kDashed))
159frame4.SetMinimum(0)
160
161c = ROOT.TCanvas("rf610_visualerror", "rf610_visualerror", 800, 800)
162c.Divide(2, 2)
163c.cd(1)
164ROOT.gPad.SetLeftMargin(0.15)
165frame.GetYaxis().SetTitleOffset(1.4)
166frame.Draw()
167c.cd(2)
168ROOT.gPad.SetLeftMargin(0.15)
169frame2.GetYaxis().SetTitleOffset(1.6)
170frame2.Draw()
171c.cd(3)
172ROOT.gPad.SetLeftMargin(0.15)
173frame3.GetYaxis().SetTitleOffset(1.6)
174frame3.Draw()
175c.cd(4)
176ROOT.gPad.SetLeftMargin(0.15)
177frame4.GetYaxis().SetTitleOffset(1.6)
178frame4.Draw()
179
180c.SaveAs("rf610_visualerror.png")