rf707_kernelestimation.cxx

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 "RooPolynomial.h"
00020 #include "RooKeysPdf.h"
00021 #include "RooNDKeysPdf.h"
00022 #include "RooProdPdf.h"
00023 #include "TCanvas.h"
00024 #include "TH1.h"
00025 #include "RooPlot.h"
00026 using namespace RooFit ;
00027 
00028 
00029 // Elementary operations on a gaussian PDF
00030 class TestBasic707 : public RooFitTestUnit
00031 {
00032 public: 
00033   TestBasic707(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Kernel estimation p.d.f.s",refFile,writeRef,verbose) {} ;
00034   Bool_t testCode() {
00035 
00036   // C r e a t e   l o w   s t a t s   1 - D   d a t a s e t 
00037   // -------------------------------------------------------
00038 
00039   // Create a toy pdf for sampling
00040   RooRealVar x("x","x",0,20) ;
00041   RooPolynomial p("p","p",x,RooArgList(RooConst(0.01),RooConst(-0.01),RooConst(0.0004))) ;
00042 
00043   // Sample 500 events from p
00044   RooDataSet* data1 = p.generate(x,200) ;
00045 
00046 
00047 
00048   // 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
00049   // ---------------------------------------------------------------
00050 
00051   // Create adaptive kernel estimation pdf. In this configuration the input data
00052   // is mirrored over the boundaries to minimize edge effects in distribution
00053   // that do not fall to zero towards the edges
00054   RooKeysPdf kest1("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth) ;
00055 
00056   // An adaptive kernel estimation pdf on the same data without mirroring option
00057   // for comparison
00058   RooKeysPdf kest2("kest2","kest2",x,*data1,RooKeysPdf::NoMirror) ;
00059 
00060 
00061   // Adaptive kernel estimation pdf with increased bandwidth scale factor
00062   // (promotes smoothness over detail preservation)
00063   RooKeysPdf kest3("kest3","kest3",x,*data1,RooKeysPdf::MirrorBoth,2) ;
00064 
00065 
00066   // Plot kernel estimation pdfs with and without mirroring over data
00067   RooPlot* frame = x.frame(Title("Adaptive kernel estimation pdf with and w/o mirroring"),Bins(20)) ;
00068   data1->plotOn(frame) ;
00069   kest1.plotOn(frame) ;    
00070   kest2.plotOn(frame,LineStyle(kDashed),LineColor(kRed)) ;    
00071 
00072 
00073   // Plot kernel estimation pdfs with regular and increased bandwidth
00074   RooPlot* frame2 = x.frame(Title("Adaptive kernel estimation pdf with regular, increased bandwidth")) ;
00075   kest1.plotOn(frame2) ;    
00076   kest3.plotOn(frame2,LineColor(kMagenta)) ;    
00077 
00078 
00079 
00080   // C r e a t e   l o w   s t a t s   2 - D   d a t a s e t 
00081   // -------------------------------------------------------
00082 
00083   // Construct a 2D toy pdf for sampleing
00084   RooRealVar y("y","y",0,20) ;
00085   RooPolynomial py("py","py",y,RooArgList(RooConst(0.01),RooConst(0.01),RooConst(-0.0004))) ;
00086   RooProdPdf pxy("pxy","pxy",RooArgSet(p,py)) ;
00087   RooDataSet* data2 = pxy.generate(RooArgSet(x,y),1000) ;
00088 
00089 
00090 
00091   // 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
00092   // ---------------------------------------------------------------
00093 
00094   // Create 2D adaptive kernel estimation pdf with mirroring 
00095   RooNDKeysPdf kest4("kest4","kest4",RooArgSet(x,y),*data2,"am") ;
00096 
00097   // Create 2D adaptive kernel estimation pdf with mirroring and double bandwidth
00098   RooNDKeysPdf kest5("kest5","kest5",RooArgSet(x,y),*data2,"am",2) ;
00099 
00100   // Create a histogram of the data
00101   TH1* hh_data = data2->createHistogram("hh_data",x,Binning(10),YVar(y,Binning(10))) ;
00102 
00103   // Create histogram of the 2d kernel estimation pdfs
00104   TH1* hh_pdf = kest4.createHistogram("hh_pdf",x,Binning(25),YVar(y,Binning(25))) ;
00105   TH1* hh_pdf2 = kest5.createHistogram("hh_pdf2",x,Binning(25),YVar(y,Binning(25))) ;
00106   hh_pdf->SetLineColor(kBlue) ;
00107   hh_pdf2->SetLineColor(kMagenta) ;
00108     
00109   regPlot(frame,"rf707_plot1") ;
00110   regPlot(frame2,"rf707_plot2") ;
00111   regTH(hh_data,"rf707_hhdata") ;
00112   regTH(hh_pdf,"rf707_hhpdf") ;
00113   regTH(hh_pdf2,"rf707_hhpdf2") ;
00114   
00115   delete data1 ;
00116   delete data2 ;
00117 
00118   return kTRUE ;
00119   }
00120 } ;

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