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_code
8/// \macro_output
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooPolynomial.h"
17#include "RooKeysPdf.h"
18#include "RooNDKeysPdf.h"
19#include "RooProdPdf.h"
20#include "TCanvas.h"
21#include "TAxis.h"
22#include "TH1.h"
23#include "RooPlot.h"
24using namespace RooFit;
25
27{
28 // C r e a t e l o w s t a t s 1 - D d a t a s e t
29 // -------------------------------------------------------
30
31 // Create a toy pdf for sampling
32 RooRealVar x("x", "x", 0, 20);
33 RooPolynomial p("p", "p", x, RooArgList(0.01, -0.01, 0.0004));
34
35 // Sample 500 events from p
36 RooDataSet *data1 = p.generate(x, 200);
37
38 // 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
39 // ---------------------------------------------------------------
40
41 // Create adaptive kernel estimation pdf. In this configuration the input data
42 // is mirrored over the boundaries to minimize edge effects in distribution
43 // that do not fall to zero towards the edges
44 RooKeysPdf kest1("kest1", "kest1", x, *data1, RooKeysPdf::MirrorBoth);
45
46 // An adaptive kernel estimation pdf on the same data without mirroring option
47 // for comparison
48 RooKeysPdf kest2("kest2", "kest2", x, *data1, RooKeysPdf::NoMirror);
49
50 // Adaptive kernel estimation pdf with increased bandwidth scale factor
51 // (promotes smoothness over detail preservation)
52 RooKeysPdf kest3("kest1", "kest1", x, *data1, RooKeysPdf::MirrorBoth, 2);
53
54 // Plot kernel estimation pdfs with and without mirroring over data
55 RooPlot *frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"), Bins(20));
56 data1->plotOn(frame);
57 kest1.plotOn(frame);
58 kest2.plotOn(frame, LineStyle(kDashed), LineColor(kRed));
59
60 // Plot kernel estimation pdfs with regular and increased bandwidth
61 RooPlot *frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth"));
62 kest1.plotOn(frame2);
63 kest3.plotOn(frame2, LineColor(kMagenta));
64
65 // C r e a t e l o w s t a t s 2 - D d a t a s e t
66 // -------------------------------------------------------
67
68 // Construct a 2D toy pdf for sampling
69 RooRealVar y("y", "y", 0, 20);
70 RooPolynomial py("py", "py", y, RooArgList(0.01, 0.01, -0.0004));
71 RooProdPdf pxy("pxy", "pxy", RooArgSet(p, py));
72 RooDataSet *data2 = pxy.generate(RooArgSet(x, y), 1000);
73
74 // 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
75 // ---------------------------------------------------------------
76
77 // Create 2D adaptive kernel estimation pdf with mirroring
78 RooNDKeysPdf kest4("kest4", "kest4", RooArgSet(x, y), *data2, "am");
79
80 // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
81 RooNDKeysPdf kest5("kest5", "kest5", RooArgSet(x, y), *data2, "am", 2);
82
83 // Create a histogram of the data
84 TH1 *hh_data = data2->createHistogram("hh_data", x, Binning(10), YVar(y, Binning(10)));
85
86 // Create histogram of the 2d kernel estimation pdfs
87 TH1 *hh_pdf = kest4.createHistogram("hh_pdf", x, Binning(25), YVar(y, Binning(25)));
88 TH1 *hh_pdf2 = kest5.createHistogram("hh_pdf2", x, Binning(25), YVar(y, Binning(25)));
89 hh_pdf->SetLineColor(kBlue);
90 hh_pdf2->SetLineColor(kMagenta);
91
92 TCanvas *c = new TCanvas("rf707_kernelestimation", "rf707_kernelestimation", 800, 800);
93 c->Divide(2, 2);
94 c->cd(1);
95 gPad->SetLeftMargin(0.15);
96 frame->GetYaxis()->SetTitleOffset(1.4);
97 frame->Draw();
98 c->cd(2);
99 gPad->SetLeftMargin(0.15);
100 frame2->GetYaxis()->SetTitleOffset(1.8);
101 frame2->Draw();
102 c->cd(3);
103 gPad->SetLeftMargin(0.15);
104 hh_data->GetZaxis()->SetTitleOffset(1.4);
105 hh_data->Draw("lego");
106 c->cd(4);
107 gPad->SetLeftMargin(0.20);
108 hh_pdf->GetZaxis()->SetTitleOffset(2.4);
109 hh_pdf->Draw("surf");
110 hh_pdf2->Draw("surfsame");
111}
#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
winID h TVirtualViewer3D TVirtualGLPainter p
#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
TH1 * createHistogram(const char *name, const RooAbsRealLValue &xvar, 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
Calls createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooLinkedList& argList) c...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
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:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:239
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
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:33
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
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:324
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
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...
Definition Common.h:18
const char * Title
Definition TXMLSetup.cxx:68