ROOT logo

From $ROOTSYS/tutorials/roofit/rf707_kernelestimation.C

//////////////////////////////////////////////////////////////////////////
//
// 'SPECIAL PDFS' RooFit tutorial macro #707
// 
// Using non-parametric (multi-dimensional) kernel estimation p.d.f.s
//
//
//
// 07/2008 - Wouter Verkerke 
// 
/////////////////////////////////////////////////////////////////////////

#ifndef __CINT__
#include "RooGlobalFunc.h"
#endif
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooPolynomial.h"
#include "RooKeysPdf.h"
#include "RooNDKeysPdf.h"
#include "RooProdPdf.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
#include "RooPlot.h"
using namespace RooFit ;


void rf707_kernelestimation()
{
  // C r e a t e   l o w   s t a t s   1 - D   d a t a s e t 
  // -------------------------------------------------------

  // Create a toy pdf for sampling
  RooRealVar x("x","x",0,20) ;
  RooPolynomial p("p","p",x,RooArgList(RooConst(0.01),RooConst(-0.01),RooConst(0.0004))) ;

  // Sample 500 events from p
  RooDataSet* data1 = p.generate(x,200) ;



  // 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
  // ---------------------------------------------------------------

  // Create adaptive kernel estimation pdf. In this configuration the input data
  // is mirrored over the boundaries to minimize edge effects in distribution
  // that do not fall to zero towards the edges
  RooKeysPdf kest1("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth) ;

  // An adaptive kernel estimation pdf on the same data without mirroring option
  // for comparison
  RooKeysPdf kest2("kest2","kest2",x,*data1,RooKeysPdf::NoMirror) ;


  // Adaptive kernel estimation pdf with increased bandwidth scale factor
  // (promotes smoothness over detail preservation)
  RooKeysPdf kest3("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth,2) ;


  // Plot kernel estimation pdfs with and without mirroring over data
  RooPlot* frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"),Bins(20)) ;
  data1->plotOn(frame) ;
  kest1.plotOn(frame) ;    
  kest2.plotOn(frame,LineStyle(kDashed),LineColor(kRed)) ;    


  // Plot kernel estimation pdfs with regular and increased bandwidth
  RooPlot* frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth")) ;
  kest1.plotOn(frame2) ;    
  kest3.plotOn(frame2,LineColor(kMagenta)) ;    



  // C r e a t e   l o w   s t a t s   2 - D   d a t a s e t 
  // -------------------------------------------------------

  // Construct a 2D toy pdf for sampleing
  RooRealVar y("y","y",0,20) ;
  RooPolynomial py("py","py",y,RooArgList(RooConst(0.01),RooConst(0.01),RooConst(-0.0004))) ;
  RooProdPdf pxy("pxy","pxy",RooArgSet(p,py)) ;
  RooDataSet* data2 = pxy.generate(RooArgSet(x,y),1000) ;



  // 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
  // ---------------------------------------------------------------

  // Create 2D adaptive kernel estimation pdf with mirroring 
  RooNDKeysPdf kest4("kest4","kest4",RooArgSet(x,y),*data2,"am") ;

  // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
  RooNDKeysPdf kest5("kest5","kest5",RooArgSet(x,y),*data2,"am",2) ;

  // Create a histogram of the data
  TH1* hh_data = data2->createHistogram("hh_data",x,Binning(10),YVar(y,Binning(10))) ;

  // Create histogram of the 2d kernel estimation pdfs
  TH1* hh_pdf = kest4.createHistogram("hh_pdf",x,Binning(25),YVar(y,Binning(25))) ;
  TH1* hh_pdf2 = kest5.createHistogram("hh_pdf2",x,Binning(25),YVar(y,Binning(25))) ;
  hh_pdf->SetLineColor(kBlue) ;
  hh_pdf2->SetLineColor(kMagenta) ;
  


  TCanvas* c = new TCanvas("rf707_kernelestimation","rf707_kernelestimation",800,800) ;
  c->Divide(2,2) ;
  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.8) ; frame2->Draw() ;
  c->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_data->GetZaxis()->SetTitleOffset(1.4) ; hh_data->Draw("lego") ;
  c->cd(4) ; gPad->SetLeftMargin(0.20) ; hh_pdf->GetZaxis()->SetTitleOffset(2.4) ; hh_pdf->Draw("surf") ; hh_pdf2->Draw("surfsame") ;
	        

}
 rf707_kernelestimation.C:1
 rf707_kernelestimation.C:2
 rf707_kernelestimation.C:3
 rf707_kernelestimation.C:4
 rf707_kernelestimation.C:5
 rf707_kernelestimation.C:6
 rf707_kernelestimation.C:7
 rf707_kernelestimation.C:8
 rf707_kernelestimation.C:9
 rf707_kernelestimation.C:10
 rf707_kernelestimation.C:11
 rf707_kernelestimation.C:12
 rf707_kernelestimation.C:13
 rf707_kernelestimation.C:14
 rf707_kernelestimation.C:15
 rf707_kernelestimation.C:16
 rf707_kernelestimation.C:17
 rf707_kernelestimation.C:18
 rf707_kernelestimation.C:19
 rf707_kernelestimation.C:20
 rf707_kernelestimation.C:21
 rf707_kernelestimation.C:22
 rf707_kernelestimation.C:23
 rf707_kernelestimation.C:24
 rf707_kernelestimation.C:25
 rf707_kernelestimation.C:26
 rf707_kernelestimation.C:27
 rf707_kernelestimation.C:28
 rf707_kernelestimation.C:29
 rf707_kernelestimation.C:30
 rf707_kernelestimation.C:31
 rf707_kernelestimation.C:32
 rf707_kernelestimation.C:33
 rf707_kernelestimation.C:34
 rf707_kernelestimation.C:35
 rf707_kernelestimation.C:36
 rf707_kernelestimation.C:37
 rf707_kernelestimation.C:38
 rf707_kernelestimation.C:39
 rf707_kernelestimation.C:40
 rf707_kernelestimation.C:41
 rf707_kernelestimation.C:42
 rf707_kernelestimation.C:43
 rf707_kernelestimation.C:44
 rf707_kernelestimation.C:45
 rf707_kernelestimation.C:46
 rf707_kernelestimation.C:47
 rf707_kernelestimation.C:48
 rf707_kernelestimation.C:49
 rf707_kernelestimation.C:50
 rf707_kernelestimation.C:51
 rf707_kernelestimation.C:52
 rf707_kernelestimation.C:53
 rf707_kernelestimation.C:54
 rf707_kernelestimation.C:55
 rf707_kernelestimation.C:56
 rf707_kernelestimation.C:57
 rf707_kernelestimation.C:58
 rf707_kernelestimation.C:59
 rf707_kernelestimation.C:60
 rf707_kernelestimation.C:61
 rf707_kernelestimation.C:62
 rf707_kernelestimation.C:63
 rf707_kernelestimation.C:64
 rf707_kernelestimation.C:65
 rf707_kernelestimation.C:66
 rf707_kernelestimation.C:67
 rf707_kernelestimation.C:68
 rf707_kernelestimation.C:69
 rf707_kernelestimation.C:70
 rf707_kernelestimation.C:71
 rf707_kernelestimation.C:72
 rf707_kernelestimation.C:73
 rf707_kernelestimation.C:74
 rf707_kernelestimation.C:75
 rf707_kernelestimation.C:76
 rf707_kernelestimation.C:77
 rf707_kernelestimation.C:78
 rf707_kernelestimation.C:79
 rf707_kernelestimation.C:80
 rf707_kernelestimation.C:81
 rf707_kernelestimation.C:82
 rf707_kernelestimation.C:83
 rf707_kernelestimation.C:84
 rf707_kernelestimation.C:85
 rf707_kernelestimation.C:86
 rf707_kernelestimation.C:87
 rf707_kernelestimation.C:88
 rf707_kernelestimation.C:89
 rf707_kernelestimation.C:90
 rf707_kernelestimation.C:91
 rf707_kernelestimation.C:92
 rf707_kernelestimation.C:93
 rf707_kernelestimation.C:94
 rf707_kernelestimation.C:95
 rf707_kernelestimation.C:96
 rf707_kernelestimation.C:97
 rf707_kernelestimation.C:98
 rf707_kernelestimation.C:99
 rf707_kernelestimation.C:100
 rf707_kernelestimation.C:101
 rf707_kernelestimation.C:102
 rf707_kernelestimation.C:103
 rf707_kernelestimation.C:104
 rf707_kernelestimation.C:105
 rf707_kernelestimation.C:106
 rf707_kernelestimation.C:107
 rf707_kernelestimation.C:108
 rf707_kernelestimation.C:109
 rf707_kernelestimation.C:110
 rf707_kernelestimation.C:111
 rf707_kernelestimation.C:112
 rf707_kernelestimation.C:113
 rf707_kernelestimation.C:114
 rf707_kernelestimation.C:115
 rf707_kernelestimation.C:116
 rf707_kernelestimation.C:117