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