ROOT
tags/v6-34-10
Reference Guide
Loading...
Searching...
No Matches
rf403_weightedevts.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## 'DATA AND CATEGORIES' RooFit tutorial macro #403
5
##
6
## Using weights in unbinned datasets
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C version)
14
15
from
__future__
import
print_function
16
import
ROOT
17
18
19
# Create observable and unweighted dataset
20
# -------------------------------------------
21
22
# Declare observable
23
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
24
x.setBins
(40)
25
26
# Construction a uniform pdf
27
p0 =
ROOT.RooPolynomial
(
"px"
,
"px"
, x)
28
29
# Sample 1000 events from pdf
30
data =
p0.generate
({x}, 1000)
31
32
# Calculate weight and make dataset weighted
33
# --------------------------------------------------
34
35
# Construct formula to calculate (fake) weight for events
36
wFunc =
ROOT.RooFormulaVar
(
"w"
,
"event weight"
,
"(x*x+10)"
, [x])
37
38
# Add column with variable w to previously generated dataset
39
w =
data.addColumn
(wFunc)
40
41
# Dataset d is now a dataset with two observable (x,w) with 1000 entries
42
data.Print
()
43
44
# Instruct dataset wdata in interpret w as event weight rather than as
45
# observable
46
wdata =
ROOT.RooDataSet
(
data.GetName
(),
data.GetTitle
(), data,
data.get
(),
""
,
w.GetName
())
47
48
# Dataset d is now a dataset with one observable (x) with 1000 entries and
49
# a sum of weights of ~430K
50
wdata.Print
()
51
52
# Unbinned ML fit to weighted data
53
# ---------------------------------------------------------------
54
55
# Construction quadratic polynomial pdf for fitting
56
a0 =
ROOT.RooRealVar
(
"a0"
,
"a0"
, 1)
57
a1 =
ROOT.RooRealVar
(
"a1"
,
"a1"
, 0, -1, 1)
58
a2 =
ROOT.RooRealVar
(
"a2"
,
"a2"
, 1, 0, 10)
59
p2 =
ROOT.RooPolynomial
(
"p2"
,
"p2"
, x, [a0, a1, a2], 0)
60
61
# Fit quadratic polynomial to weighted data
62
63
# NOTE: A plain Maximum likelihood fit to weighted data does in general
64
# NOT result in correct error estimates, individual
65
# event weights represent Poisson statistics themselves.
66
#
67
# Fit with 'wrong' errors
68
r_ml_wgt =
p2.fitTo
(wdata, Save=
True
, PrintLevel=-1)
69
70
# A first order correction to estimated parameter errors in an
71
# (unbinned) ML fit can be obtained by calculating the
72
# covariance matrix as
73
#
74
# V' = V C-1 V
75
#
76
# where V is the covariance matrix calculated from a fit
77
# to -logL = - sum [ w_i log f(x_i) ] and C is the covariance
78
# matrix calculated from -logL' = -sum [ w_i^2 log f(x_i) ]
79
# (i.e. the weights are applied squared)
80
#
81
# A fit in self mode can be performed as follows:
82
83
r_ml_wgt_corr =
p2.fitTo
(wdata, Save=
True
, SumW2Error=
True
, PrintLevel=-1)
84
85
# Plot weighted data and fit result
86
# ---------------------------------------------------------------
87
88
# Construct plot frame
89
frame =
x.frame
(Title=
"Unbinned ML fit, chi^2 fit to weighted data"
)
90
91
# Plot data using sum-of-weights-squared error rather than Poisson errors
92
wdata.plotOn
(frame, DataError=
"SumW2"
)
93
94
# Overlay result of 2nd order polynomial fit to weighted data
95
p2.plotOn
(frame)
96
97
# ML fit of pdf to equivalent unweighted dataset
98
# ---------------------------------------------------------------------
99
100
# Construct a pdf with the same shape as p0 after weighting
101
genPdf =
ROOT.RooGenericPdf
(
"genPdf"
,
"x*x+10"
, [x])
102
103
# Sample a dataset with the same number of events as data
104
data2 =
genPdf.generate
({x}, 1000)
105
106
# Sample a dataset with the same number of weights as data
107
data3 =
genPdf.generate
({x}, 43000)
108
109
# Fit the 2nd order polynomial to both unweighted datasets and save the
110
# results for comparison
111
r_ml_unw10 =
p2.fitTo
(data2, Save=
True
, PrintLevel=-1)
112
r_ml_unw43 =
p2.fitTo
(data3, Save=
True
, PrintLevel=-1)
113
114
# Chis2 fit of pdf to binned weighted dataset
115
# ---------------------------------------------------------------------------
116
117
# Construct binned clone of unbinned weighted dataset
118
binnedData =
wdata.binnedClone
()
119
binnedData.Print
(
"v"
)
120
121
# Perform chi2 fit to binned weighted dataset using sum-of-weights errors
122
#
123
# NB: Within the usual approximations of a chi2 fit, chi2 fit to weighted
124
# data using sum-of-weights-squared errors does give correct error
125
# estimates
126
chi2 =
p2.createChi2
(binnedData,
ROOT.RooFit.DataError
(
"SumW2"
))
127
m =
ROOT.RooMinimizer
(chi2)
128
m.migrad
()
129
m.hesse
()
130
131
# Plot chi^2 fit result on frame as well
132
r_chi2_wgt =
m.save
()
133
p2.plotOn
(frame, LineStyle=
"--"
, LineColor=
"r"
)
134
135
# Compare fit results of chi2, L fits to (un)weighted data
136
# ------------------------------------------------------------
137
138
# Note that ML fit on 1Kevt of weighted data is closer to result of ML fit on 43Kevt of unweighted data
139
# than to 1Kevt of unweighted data, the reference chi^2 fit with SumW2 error gives a result closer to
140
# that of an unbinned ML fit to 1Kevt of unweighted data.
141
142
print(
"==> ML Fit results on 1K unweighted events"
)
143
r_ml_unw10.Print
()
144
print(
"==> ML Fit results on 43K unweighted events"
)
145
r_ml_unw43.Print
()
146
print(
"==> ML Fit results on 1K weighted events with a summed weight of 43K"
)
147
r_ml_wgt.Print
()
148
print(
"==> Corrected ML Fit results on 1K weighted events with a summed weight of 43K"
)
149
r_ml_wgt_corr.Print
()
150
print(
"==> Chi2 Fit results on 1K weighted events with a summed weight of 43K"
)
151
r_chi2_wgt.Print
()
152
153
c =
ROOT.TCanvas
(
"rf403_weightedevts"
,
"rf403_weightedevts"
, 600, 600)
154
ROOT.gPad.SetLeftMargin
(0.15)
155
frame.GetYaxis
().SetTitleOffset(1.8)
156
frame.Draw
()
157
158
c.SaveAs
(
"rf403_weightedevts.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
rf403_weightedevts.py
ROOT tags/v6-34-10 - Reference Guide Generated on Mon Jun 30 2025 10:56:37 (GVA Time) using Doxygen 1.10.0