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