ROOT
Version v6.32
master
v6.34
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
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_image
8
## \macro_code
9
## \macro_output
10
##
11
## \date February 2018
12
## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14
import
ROOT
15
16
17
# Create model and dataset
18
# -----------------------------------------------
19
20
# Observable
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -20, 20)
22
23
# Model (intentional strong correlations)
24
mean =
ROOT.RooRealVar
(
"mean"
,
"mean of g1 and g2"
, 0, -1, 1)
25
sigma_g1 =
ROOT.RooRealVar
(
"sigma_g1"
,
"width of g1"
, 2)
26
g1 =
ROOT.RooGaussian
(
"g1"
,
"g1"
, x, mean, sigma_g1)
27
28
sigma_g2 =
ROOT.RooRealVar
(
"sigma_g2"
,
"width of g2"
, 4, 3.0, 5.0)
29
g2 =
ROOT.RooGaussian
(
"g2"
,
"g2"
, x, mean, sigma_g2)
30
31
frac =
ROOT.RooRealVar
(
"frac"
,
"frac"
, 0.5, 0.0, 1.0)
32
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [g1, g2], [frac])
33
34
# Generate 1000 events
35
data =
model.generate
({x}, 1000)
36
37
# Fit model to data
38
# ----------------------------------
39
40
r =
model.fitTo
(data, Save=
True
, PrintLevel=-1)
41
42
# Create MV Gaussian pdf of fitted parameters
43
# ------------------------------------------------------------------------------------
44
45
parabPdf =
r.createHessePdf
({frac, mean, sigma_g2})
46
47
# Some exercises with the parameter pdf
48
# -----------------------------------------------------------------------------
49
50
# Generate 100K points in the parameter space, from the MVGaussian pdf
51
d =
parabPdf.generate
({mean, sigma_g2, frac}, 100000)
52
53
# Sample a 3-D histogram of the pdf to be visualized as an error
54
# ellipsoid using the GLISO draw option
55
hh_3d =
parabPdf.createHistogram
(
"mean,sigma_g2,frac"
, 25, 25, 25)
56
hh_3d.SetFillColor
(
ROOT.kBlue
)
57
58
# Project 3D parameter pdf down to 3 permutations of two-dimensional pdfs
59
# The integrations corresponding to these projections are performed analytically
60
# by the MV Gaussian pdf
61
pdf_sigmag2_frac =
parabPdf.createProjection
({mean})
62
pdf_mean_frac =
parabPdf.createProjection
({sigma_g2})
63
pdf_mean_sigmag2 =
parabPdf.createProjection
({frac})
64
65
# Make 2D plots of the 3 two-dimensional pdf projections
66
hh_sigmag2_frac =
pdf_sigmag2_frac.createHistogram
(
"sigma_g2,frac"
, 50, 50)
67
hh_mean_frac =
pdf_mean_frac.createHistogram
(
"mean,frac"
, 50, 50)
68
hh_mean_sigmag2 =
pdf_mean_sigmag2.createHistogram
(
"mean,sigma_g2"
, 50, 50)
69
hh_mean_frac.SetLineColor
(
ROOT.kBlue
)
70
hh_sigmag2_frac.SetLineColor
(
ROOT.kBlue
)
71
hh_mean_sigmag2.SetLineColor
(
ROOT.kBlue
)
72
73
# Draw the 'sigar'
74
ROOT.gStyle.SetCanvasPreferGL
(
True
)
75
ROOT.gStyle.SetPalette
(1)
76
c1 =
ROOT.TCanvas
(
"rf608_fitresultaspdf_1"
,
"rf608_fitresultaspdf_1"
, 600, 600)
77
hh_3d.Draw
(
"gliso"
)
78
79
c1.SaveAs
(
"rf608_fitresultaspdf_1.png"
)
80
81
# Draw the 2D projections of the 3D pdf
82
c2 =
ROOT.TCanvas
(
"rf608_fitresultaspdf_2"
,
"rf608_fitresultaspdf_2"
, 900, 600)
83
c2.Divide
(3, 2)
84
c2.cd
(1)
85
ROOT.gPad.SetLeftMargin
(0.15)
86
hh_mean_sigmag2.GetZaxis
().SetTitleOffset(1.4)
87
hh_mean_sigmag2.Draw
(
"surf3"
)
88
c2.cd
(2)
89
ROOT.gPad.SetLeftMargin
(0.15)
90
hh_sigmag2_frac.GetZaxis
().SetTitleOffset(1.4)
91
hh_sigmag2_frac.Draw
(
"surf3"
)
92
c2.cd
(3)
93
ROOT.gPad.SetLeftMargin
(0.15)
94
hh_mean_frac.GetZaxis
().SetTitleOffset(1.4)
95
hh_mean_frac.Draw
(
"surf3"
)
96
97
# Draw the distributions of parameter points sampled from the pdf
98
tmp1 =
d.createHistogram
(mean, sigma_g2, 50, 50)
99
tmp2 =
d.createHistogram
(sigma_g2, frac, 50, 50)
100
tmp3 =
d.createHistogram
(mean, frac, 50, 50)
101
102
c2.cd
(4)
103
ROOT.gPad.SetLeftMargin
(0.15)
104
tmp1.GetZaxis
().SetTitleOffset(1.4)
105
tmp1.Draw
(
"lego3"
)
106
c2.cd
(5)
107
ROOT.gPad.SetLeftMargin
(0.15)
108
tmp2.GetZaxis
().SetTitleOffset(1.4)
109
tmp2.Draw
(
"lego3"
)
110
c2.cd
(6)
111
ROOT.gPad.SetLeftMargin
(0.15)
112
tmp3.GetZaxis
().SetTitleOffset(1.4)
113
tmp3.Draw
(
"lego3"
)
114
115
c2.SaveAs
(
"rf608_fitresultaspdf_2.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
rf608_fitresultaspdf.py
ROOT v6-32 - Reference Guide Generated on Wed Apr 2 2025 08:24:49 (GVA Time) using Doxygen 1.10.0