Logo ROOT   6.16/01
Reference Guide
rf707_kernelestimation.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'SPECIAL PDFS' RooFit tutorial macro #707
6##
7## Using non-parametric (multi-dimensional) kernel estimation p.d.f.s
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13
14
15import ROOT
16
17
18# Create low stats 1D dataset
19# -------------------------------------------------------
20
21# Create a toy pdf for sampling
22x = ROOT.RooRealVar("x", "x", 0, 20)
23p = ROOT.RooPolynomial("p", "p", x, ROOT.RooArgList(ROOT.RooFit.RooConst(
24 0.01), ROOT.RooFit.RooConst(-0.01), ROOT.RooFit.RooConst(0.0004)))
25
26# Sample 500 events from p
27data1 = p.generate(ROOT.RooArgSet(x), 200)
28
29# Create 1D kernel estimation pdf
30# ---------------------------------------------------------------
31
32# Create adaptive kernel estimation pdf. In self configuration the input data
33# is mirrored over the boundaries to minimize edge effects in distribution
34# that do not fall to zero towards the edges
35kest1 = ROOT.RooKeysPdf("kest1", "kest1", x, data1,
36 ROOT.RooKeysPdf.MirrorBoth)
37
38# An adaptive kernel estimation pdf on the same data without mirroring option
39# for comparison
40kest2 = ROOT.RooKeysPdf("kest2", "kest2", x, data1,
41 ROOT.RooKeysPdf.NoMirror)
42
43# Adaptive kernel estimation pdf with increased bandwidth scale factor
44# (promotes smoothness over detail preservation)
45kest3 = ROOT.RooKeysPdf("kest1", "kest1", x, data1,
46 ROOT.RooKeysPdf.MirrorBoth, 2)
47
48# Plot kernel estimation pdfs with and without mirroring over data
49frame = x.frame(ROOT.RooFit.Title(
50 "Adaptive kernel estimation pdf with and w/o mirroring"), ROOT.RooFit.Bins(20))
51data1.plotOn(frame)
52kest1.plotOn(frame)
53kest2.plotOn(frame, ROOT.RooFit.LineStyle(
54 ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
55
56# Plot kernel estimation pdfs with regular and increased bandwidth
57frame2 = x.frame(ROOT.RooFit.Title(
58 "Adaptive kernel estimation pdf with regular, bandwidth"))
59kest1.plotOn(frame2)
60kest3.plotOn(frame2, ROOT.RooFit.LineColor(ROOT.kMagenta))
61
62# Create low status 2D dataset
63# -------------------------------------------------------
64
65# Construct a 2D toy pdf for sampleing
66y = ROOT.RooRealVar("y", "y", 0, 20)
67py = ROOT.RooPolynomial("py", "py", y, ROOT.RooArgList(ROOT.RooFit.RooConst(0.01), ROOT.RooFit.RooConst(0.01), ROOT.RooFit.RooConst(-0.0004)))
68pxy = ROOT.RooProdPdf("pxy", "pxy", ROOT.RooArgList(p, py))
69data2 = pxy.generate(ROOT.RooArgSet(x, y), 1000)
70
71# Create 2D kernel estimation pdf
72# ---------------------------------------------------------------
73
74# Create 2D adaptive kernel estimation pdf with mirroring
75kest4 = ROOT.RooNDKeysPdf("kest4", "kest4", ROOT.RooArgList(x, y), data2, "am")
76
77# Create 2D adaptive kernel estimation pdf with mirroring and double
78# bandwidth
79kest5 = ROOT.RooNDKeysPdf("kest5", "kest5", ROOT.RooArgList(x, y), data2, "am", 2)
80
81# Create a histogram of the data
82hh_data = ROOT.RooAbsData.createHistogram(data2, "hh_data", x, ROOT.RooFit.Binning(
83 10), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(10)))
84
85# Create histogram of the 2d kernel estimation pdfs
86hh_pdf = kest4.createHistogram("hh_pdf", x, ROOT.RooFit.Binning(
87 25), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(25)))
88hh_pdf2 = kest5.createHistogram("hh_pdf2", x, ROOT.RooFit.Binning(
89 25), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(25)))
90hh_pdf.SetLineColor(ROOT.kBlue)
91hh_pdf2.SetLineColor(ROOT.kMagenta)
92
93c = ROOT.TCanvas("rf707_kernelestimation",
94 "rf707_kernelestimation", 800, 800)
95c.Divide(2, 2)
96c.cd(1)
97ROOT.gPad.SetLeftMargin(0.15)
98frame.GetYaxis().SetTitleOffset(1.4)
99frame.Draw()
100c.cd(2)
101ROOT.gPad.SetLeftMargin(0.15)
102frame2.GetYaxis().SetTitleOffset(1.8)
103frame2.Draw()
104c.cd(3)
105ROOT.gPad.SetLeftMargin(0.15)
106hh_data.GetZaxis().SetTitleOffset(1.4)
107hh_data.Draw("lego")
108c.cd(4)
109ROOT.gPad.SetLeftMargin(0.20)
110hh_pdf.GetZaxis().SetTitleOffset(2.4)
111hh_pdf.Draw("surf")
112hh_pdf2.Draw("surfsame")
113
114c.SaveAs("rf707_kernelestimation.png")