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