Logo ROOT   6.16/01
Reference Guide
rf608_fitresultaspdf.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #608
6##
7## Representing the parabolic approximation of the fit as
8# a multi-variate Gaussian on the parameters of the fitted p.d.f.
9#
10##
11## \macro_code
12##
13## \date February 2018
14## \author Clemens Lange
15## \author Wouter Verkerke (C version)
16
17
18import ROOT
19
20
21# Create model and dataset
22# -----------------------------------------------
23
24# Observable
25x = ROOT.RooRealVar("x", "x", -20, 20)
26
27# Model (intentional strong correlations)
28mean = ROOT.RooRealVar("mean", "mean of g1 and g2", 0, -1, 1)
29sigma_g1 = ROOT.RooRealVar("sigma_g1", "width of g1", 2)
30g1 = ROOT.RooGaussian("g1", "g1", x, mean, sigma_g1)
31
32sigma_g2 = ROOT.RooRealVar("sigma_g2", "width of g2", 4, 3.0, 5.0)
33g2 = ROOT.RooGaussian("g2", "g2", x, mean, sigma_g2)
34
35frac = ROOT.RooRealVar("frac", "frac", 0.5, 0.0, 1.0)
36model = ROOT.RooAddPdf("model", "model", ROOT.RooArgList(g1, g2), ROOT.RooArgList(frac))
37
38# Generate 1000 events
39data = model.generate(ROOT.RooArgSet(x), 1000)
40
41# Fit model to data
42# ----------------------------------
43
44r = model.fitTo(data, ROOT.RooFit.Save())
45
46# Create MV Gaussian pdf of fitted parameters
47# ------------------------------------------------------------------------------------
48
49parabPdf = r.createHessePdf(ROOT.RooArgSet(frac, mean, sigma_g2))
50
51# Some exercises with the parameter pdf
52# -----------------------------------------------------------------------------
53
54# Generate 100K points in the parameter space, from the MVGaussian p.d.f.
55d = parabPdf.generate(ROOT.RooArgSet(mean, sigma_g2, frac), 100000)
56
57# Sample a 3-D histogram of the p.d.f. to be visualized as an error
58# ellipsoid using the GLISO draw option
59hh_3d = parabPdf.createHistogram("mean,sigma_g2,frac", 25, 25, 25)
60hh_3d.SetFillColor(ROOT.kBlue)
61
62# Project 3D parameter p.d.f. down to 3 permutations of two-dimensional p.d.f.s
63# The integrations corresponding to these projections are performed analytically
64# by the MV Gaussian p.d.f.
65pdf_sigmag2_frac = parabPdf.createProjection(ROOT.RooArgSet(mean))
66pdf_mean_frac = parabPdf.createProjection(ROOT.RooArgSet(sigma_g2))
67pdf_mean_sigmag2 = parabPdf.createProjection(ROOT.RooArgSet(frac))
68
69# Make 2D plots of the 3 two-dimensional p.d.f. projections
70hh_sigmag2_frac = pdf_sigmag2_frac.createHistogram("sigma_g2,frac", 50, 50)
71hh_mean_frac = pdf_mean_frac.createHistogram("mean,frac", 50, 50)
72hh_mean_sigmag2 = pdf_mean_sigmag2.createHistogram("mean,sigma_g2", 50, 50)
73hh_mean_frac.SetLineColor(ROOT.kBlue)
74hh_sigmag2_frac.SetLineColor(ROOT.kBlue)
75hh_mean_sigmag2.SetLineColor(ROOT.kBlue)
76
77# Draw the 'sigar'
78ROOT.gStyle.SetCanvasPreferGL(True)
79ROOT.gStyle.SetPalette(1)
80c1 = ROOT.TCanvas("rf608_fitresultaspdf_1", "rf608_fitresultaspdf_1", 600, 600)
81hh_3d.Draw("gliso")
82
83c1.SaveAs("rf608_fitresultaspdf_1.png")
84
85# Draw the 2D projections of the 3D p.d.f.
86c2 = ROOT.TCanvas("rf608_fitresultaspdf_2",
87 "rf608_fitresultaspdf_2", 900, 600)
88c2.Divide(3, 2)
89c2.cd(1)
90ROOT.gPad.SetLeftMargin(0.15)
91hh_mean_sigmag2.GetZaxis().SetTitleOffset(1.4)
92hh_mean_sigmag2.Draw("surf3")
93c2.cd(2)
94ROOT.gPad.SetLeftMargin(0.15)
95hh_sigmag2_frac.GetZaxis().SetTitleOffset(1.4)
96hh_sigmag2_frac.Draw("surf3")
97c2.cd(3)
98ROOT.gPad.SetLeftMargin(0.15)
99hh_mean_frac.GetZaxis().SetTitleOffset(1.4)
100hh_mean_frac.Draw("surf3")
101
102# Draw the distributions of parameter points sampled from the p.d.f.
103tmp1 = d.createHistogram(mean, sigma_g2, 50, 50)
104tmp2 = d.createHistogram(sigma_g2, frac, 50, 50)
105tmp3 = d.createHistogram(mean, frac, 50, 50)
106
107c2.cd(4)
108ROOT.gPad.SetLeftMargin(0.15)
109tmp1.GetZaxis().SetTitleOffset(1.4)
110tmp1.Draw("lego3")
111c2.cd(5)
112ROOT.gPad.SetLeftMargin(0.15)
113tmp2.GetZaxis().SetTitleOffset(1.4)
114tmp2.Draw("lego3")
115c2.cd(6)
116ROOT.gPad.SetLeftMargin(0.15)
117tmp3.GetZaxis().SetTitleOffset(1.4)
118tmp3.Draw("lego3")
119
120c2.SaveAs("rf608_fitresultaspdf_2.png")