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