rf707_kernelestimation.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'SPECIAL PDFS' RooFit tutorial macro #707
00004 // 
00005 // Using non-parametric (multi-dimensional) kernel estimation p.d.f.s
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussian.h"
00019 #include "RooConstVar.h"
00020 #include "RooPolynomial.h"
00021 #include "RooKeysPdf.h"
00022 #include "RooNDKeysPdf.h"
00023 #include "RooProdPdf.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "TH1.h"
00027 #include "RooPlot.h"
00028 using namespace RooFit ;
00029 
00030 
00031 void rf707_kernelestimation()
00032 {
00033   // C r e a t e   l o w   s t a t s   1 - D   d a t a s e t 
00034   // -------------------------------------------------------
00035 
00036   // Create a toy pdf for sampling
00037   RooRealVar x("x","x",0,20) ;
00038   RooPolynomial p("p","p",x,RooArgList(RooConst(0.01),RooConst(-0.01),RooConst(0.0004))) ;
00039 
00040   // Sample 500 events from p
00041   RooDataSet* data1 = p.generate(x,200) ;
00042 
00043 
00044 
00045   // 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
00046   // ---------------------------------------------------------------
00047 
00048   // Create adaptive kernel estimation pdf. In this configuration the input data
00049   // is mirrored over the boundaries to minimize edge effects in distribution
00050   // that do not fall to zero towards the edges
00051   RooKeysPdf kest1("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth) ;
00052 
00053   // An adaptive kernel estimation pdf on the same data without mirroring option
00054   // for comparison
00055   RooKeysPdf kest2("kest2","kest2",x,*data1,RooKeysPdf::NoMirror) ;
00056 
00057 
00058   // Adaptive kernel estimation pdf with increased bandwidth scale factor
00059   // (promotes smoothness over detail preservation)
00060   RooKeysPdf kest3("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth,2) ;
00061 
00062 
00063   // Plot kernel estimation pdfs with and without mirroring over data
00064   RooPlot* frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"),Bins(20)) ;
00065   data1->plotOn(frame) ;
00066   kest1.plotOn(frame) ;    
00067   kest2.plotOn(frame,LineStyle(kDashed),LineColor(kRed)) ;    
00068 
00069 
00070   // Plot kernel estimation pdfs with regular and increased bandwidth
00071   RooPlot* frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth")) ;
00072   kest1.plotOn(frame2) ;    
00073   kest3.plotOn(frame2,LineColor(kMagenta)) ;    
00074 
00075 
00076 
00077   // C r e a t e   l o w   s t a t s   2 - D   d a t a s e t 
00078   // -------------------------------------------------------
00079 
00080   // Construct a 2D toy pdf for sampleing
00081   RooRealVar y("y","y",0,20) ;
00082   RooPolynomial py("py","py",y,RooArgList(RooConst(0.01),RooConst(0.01),RooConst(-0.0004))) ;
00083   RooProdPdf pxy("pxy","pxy",RooArgSet(p,py)) ;
00084   RooDataSet* data2 = pxy.generate(RooArgSet(x,y),1000) ;
00085 
00086 
00087 
00088   // 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
00089   // ---------------------------------------------------------------
00090 
00091   // Create 2D adaptive kernel estimation pdf with mirroring 
00092   RooNDKeysPdf kest4("kest4","kest4",RooArgSet(x,y),*data2,"am") ;
00093 
00094   // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
00095   RooNDKeysPdf kest5("kest5","kest5",RooArgSet(x,y),*data2,"am",2) ;
00096 
00097   // Create a histogram of the data
00098   TH1* hh_data = data2->createHistogram("hh_data",x,Binning(10),YVar(y,Binning(10))) ;
00099 
00100   // Create histogram of the 2d kernel estimation pdfs
00101   TH1* hh_pdf = kest4.createHistogram("hh_pdf",x,Binning(25),YVar(y,Binning(25))) ;
00102   TH1* hh_pdf2 = kest5.createHistogram("hh_pdf2",x,Binning(25),YVar(y,Binning(25))) ;
00103   hh_pdf->SetLineColor(kBlue) ;
00104   hh_pdf2->SetLineColor(kMagenta) ;
00105   
00106 
00107 
00108   TCanvas* c = new TCanvas("rf707_kernelestimation","rf707_kernelestimation",800,800) ;
00109   c->Divide(2,2) ;
00110   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00111   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.8) ; frame2->Draw() ;
00112   c->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_data->GetZaxis()->SetTitleOffset(1.4) ; hh_data->Draw("lego") ;
00113   c->cd(4) ; gPad->SetLeftMargin(0.20) ; hh_pdf->GetZaxis()->SetTitleOffset(2.4) ; hh_pdf->Draw("surf") ; hh_pdf2->Draw("surfsame") ;
00114                 
00115 
00116 }

Generated on Tue Jul 5 15:45:08 2011 for ROOT_528-00b_version by  doxygen 1.5.1