ROOT  6.06/09
Reference Guide
rf707_kernelestimation.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'SPECIAL PDFS' RooFit tutorial macro #707
4 //
5 // Using non-parametric (multi-dimensional) kernel estimation p.d.f.s
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooPolynomial.h"
20 #include "RooKeysPdf.h"
21 #include "RooNDKeysPdf.h"
22 #include "RooProdPdf.h"
23 #include "TCanvas.h"
24 #include "TH1.h"
25 #include "RooPlot.h"
26 using namespace RooFit ;
27 
28 
29 // Elementary operations on a gaussian PDF
30 class TestBasic707 : public RooFitTestUnit
31 {
32 public:
33  TestBasic707(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Kernel estimation p.d.f.s",refFile,writeRef,verbose) {} ;
34  Bool_t testCode() {
35 
36  // C r e a t e l o w s t a t s 1 - D d a t a s e t
37  // -------------------------------------------------------
38 
39  // Create a toy pdf for sampling
40  RooRealVar x("x","x",0,20) ;
41  RooPolynomial p("p","p",x,RooArgList(RooConst(0.01),RooConst(-0.01),RooConst(0.0004))) ;
42 
43  // Sample 500 events from p
44  RooDataSet* data1 = p.generate(x,200) ;
45 
46 
47 
48  // 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
49  // ---------------------------------------------------------------
50 
51  // Create adaptive kernel estimation pdf. In this configuration the input data
52  // is mirrored over the boundaries to minimize edge effects in distribution
53  // that do not fall to zero towards the edges
54  RooKeysPdf kest1("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth) ;
55 
56  // An adaptive kernel estimation pdf on the same data without mirroring option
57  // for comparison
58  RooKeysPdf kest2("kest2","kest2",x,*data1,RooKeysPdf::NoMirror) ;
59 
60 
61  // Adaptive kernel estimation pdf with increased bandwidth scale factor
62  // (promotes smoothness over detail preservation)
63  RooKeysPdf kest3("kest3","kest3",x,*data1,RooKeysPdf::MirrorBoth,2) ;
64 
65 
66  // Plot kernel estimation pdfs with and without mirroring over data
67  RooPlot* frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"),Bins(20)) ;
68  data1->plotOn(frame) ;
69  kest1.plotOn(frame) ;
70  kest2.plotOn(frame,LineStyle(kDashed),LineColor(kRed)) ;
71 
72 
73  // Plot kernel estimation pdfs with regular and increased bandwidth
74  RooPlot* frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth")) ;
75  kest1.plotOn(frame2) ;
76  kest3.plotOn(frame2,LineColor(kMagenta)) ;
77 
78 
79 
80  // C r e a t e l o w s t a t s 2 - D d a t a s e t
81  // -------------------------------------------------------
82 
83  // Construct a 2D toy pdf for sampleing
84  RooRealVar y("y","y",0,20) ;
85  RooPolynomial py("py","py",y,RooArgList(RooConst(0.01),RooConst(0.01),RooConst(-0.0004))) ;
86  RooProdPdf pxy("pxy","pxy",RooArgSet(p,py)) ;
87  RooDataSet* data2 = pxy.generate(RooArgSet(x,y),1000) ;
88 
89 
90 
91  // 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
92  // ---------------------------------------------------------------
93 
94  // Create 2D adaptive kernel estimation pdf with mirroring
95  RooNDKeysPdf kest4("kest4","kest4",RooArgSet(x,y),*data2,"am") ;
96 
97  // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
98  RooNDKeysPdf kest5("kest5","kest5",RooArgSet(x,y),*data2,"am",2) ;
99 
100  // Create a histogram of the data
101  TH1* hh_data = data2->createHistogram("hh_data",x,Binning(10),YVar(y,Binning(10))) ;
102 
103  // Create histogram of the 2d kernel estimation pdfs
104  TH1* hh_pdf = kest4.createHistogram("hh_pdf",x,Binning(25),YVar(y,Binning(25))) ;
105  TH1* hh_pdf2 = kest5.createHistogram("hh_pdf2",x,Binning(25),YVar(y,Binning(25))) ;
106  hh_pdf->SetLineColor(kBlue) ;
107  hh_pdf2->SetLineColor(kMagenta) ;
108 
109  regPlot(frame,"rf707_plot1") ;
110  regPlot(frame2,"rf707_plot2") ;
111  regTH(hh_data,"rf707_hhdata") ;
112  regTH(hh_pdf,"rf707_hhpdf") ;
113  regTH(hh_pdf2,"rf707_hhpdf2") ;
114 
115  delete data1 ;
116  delete data2 ;
117 
118  return kTRUE ;
119  }
120 } ;
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:61
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
RooCmdArg LineStyle(Style_t style)
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:624
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
bool verbose
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
Double_t y[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:80
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.
RooConstVar & RooConst(Double_t val)
RooCmdArg Bins(Int_t nbin)
Definition: Rtypes.h:61
const Bool_t kTRUE
Definition: Rtypes.h:91
RooCmdArg Binning(const RooAbsBinning &binning)