Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf707_kernelestimation.C
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_output
8/// \macro_code
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooConstVar.h"
17#include "RooPolynomial.h"
18#include "RooKeysPdf.h"
19#include "RooNDKeysPdf.h"
20#include "RooProdPdf.h"
21#include "TCanvas.h"
22#include "TAxis.h"
23#include "TH1.h"
24#include "RooPlot.h"
25using namespace RooFit;
26
28{
29 // C r e a t e l o w s t a t s 1 - D d a t a s e t
30 // -------------------------------------------------------
31
32 // Create a toy pdf for sampling
33 RooRealVar x("x", "x", 0, 20);
34 RooPolynomial p("p", "p", x, RooArgList(RooConst(0.01), RooConst(-0.01), RooConst(0.0004)));
35
36 // Sample 500 events from p
37 RooDataSet *data1 = p.generate(x, 200);
38
39 // 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
40 // ---------------------------------------------------------------
41
42 // Create adaptive kernel estimation pdf. In this configuration the input data
43 // is mirrored over the boundaries to minimize edge effects in distribution
44 // that do not fall to zero towards the edges
45 RooKeysPdf kest1("kest1", "kest1", x, *data1, RooKeysPdf::MirrorBoth);
46
47 // An adaptive kernel estimation pdf on the same data without mirroring option
48 // for comparison
49 RooKeysPdf kest2("kest2", "kest2", x, *data1, RooKeysPdf::NoMirror);
50
51 // Adaptive kernel estimation pdf with increased bandwidth scale factor
52 // (promotes smoothness over detail preservation)
53 RooKeysPdf kest3("kest1", "kest1", x, *data1, RooKeysPdf::MirrorBoth, 2);
54
55 // Plot kernel estimation pdfs with and without mirroring over data
56 RooPlot *frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"), Bins(20));
57 data1->plotOn(frame);
58 kest1.plotOn(frame);
59 kest2.plotOn(frame, LineStyle(kDashed), LineColor(kRed));
60
61 // Plot kernel estimation pdfs with regular and increased bandwidth
62 RooPlot *frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth"));
63 kest1.plotOn(frame2);
64 kest3.plotOn(frame2, LineColor(kMagenta));
65
66 // C r e a t e l o w s t a t s 2 - D d a t a s e t
67 // -------------------------------------------------------
68
69 // Construct a 2D toy pdf for sampling
70 RooRealVar y("y", "y", 0, 20);
71 RooPolynomial py("py", "py", y, RooArgList(RooConst(0.01), RooConst(0.01), RooConst(-0.0004)));
72 RooProdPdf pxy("pxy", "pxy", RooArgSet(p, py));
73 RooDataSet *data2 = pxy.generate(RooArgSet(x, y), 1000);
74
75 // 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
76 // ---------------------------------------------------------------
77
78 // Create 2D adaptive kernel estimation pdf with mirroring
79 RooNDKeysPdf kest4("kest4", "kest4", RooArgSet(x, y), *data2, "am");
80
81 // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
82 RooNDKeysPdf kest5("kest5", "kest5", RooArgSet(x, y), *data2, "am", 2);
83
84 // Create a histogram of the data
85 TH1 *hh_data = data2->createHistogram("hh_data", x, Binning(10), YVar(y, Binning(10)));
86
87 // Create histogram of the 2d kernel estimation pdfs
88 TH1 *hh_pdf = kest4.createHistogram("hh_pdf", x, Binning(25), YVar(y, Binning(25)));
89 TH1 *hh_pdf2 = kest5.createHistogram("hh_pdf2", x, Binning(25), YVar(y, Binning(25)));
90 hh_pdf->SetLineColor(kBlue);
91 hh_pdf2->SetLineColor(kMagenta);
92
93 TCanvas *c = new TCanvas("rf707_kernelestimation", "rf707_kernelestimation", 800, 800);
94 c->Divide(2, 2);
95 c->cd(1);
96 gPad->SetLeftMargin(0.15);
97 frame->GetYaxis()->SetTitleOffset(1.4);
98 frame->Draw();
99 c->cd(2);
100 gPad->SetLeftMargin(0.15);
101 frame2->GetYaxis()->SetTitleOffset(1.8);
102 frame2->Draw();
103 c->cd(3);
104 gPad->SetLeftMargin(0.15);
105 hh_data->GetZaxis()->SetTitleOffset(1.4);
106 hh_data->Draw("lego");
107 c->cd(4);
108 gPad->SetLeftMargin(0.20);
109 hh_pdf->GetZaxis()->SetTitleOffset(2.4);
110 hh_pdf->Draw("surf");
111 hh_pdf2->Draw("surfsame");
112}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
#define gPad
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
TH2F * createHistogram(const RooAbsRealLValue &var1, const RooAbsRealLValue &var2, const char *cuts="", const char *name="hist") const
Create a TH2F histogram of the distribution of the specified variable using this dataset.
Class RooKeysPdf implements a one-dimensional kernel estimation p.d.f which model the distribution of...
Definition RooKeysPdf.h:25
Generic N-dimensional implementation of a kernel estimation p.d.f.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
TAxis * GetYaxis() const
Definition RooPlot.cxx:1263
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:249
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooPolynomial implements a polynomial p.d.f of the form.
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:37
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:293
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:322
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooConstVar & RooConst(Double_t val)
RooCmdArg Bins(Int_t nbin)
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition TXMLSetup.cxx:68