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