00001
00002
00003
00004
00005
00006
00007
00008
00009
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
00034
00035
00036
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
00041 RooDataSet* data1 = p.generate(x,200) ;
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051 RooKeysPdf kest1("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth) ;
00052
00053
00054
00055 RooKeysPdf kest2("kest2","kest2",x,*data1,RooKeysPdf::NoMirror) ;
00056
00057
00058
00059
00060 RooKeysPdf kest3("kest1","kest1",x,*data1,RooKeysPdf::MirrorBoth,2) ;
00061
00062
00063
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
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
00078
00079
00080
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
00089
00090
00091
00092 RooNDKeysPdf kest4("kest4","kest4",RooArgSet(x,y),*data2,"am") ;
00093
00094
00095 RooNDKeysPdf kest5("kest5","kest5",RooArgSet(x,y),*data2,"am",2) ;
00096
00097
00098 TH1* hh_data = data2->createHistogram("hh_data",x,Binning(10),YVar(y,Binning(10))) ;
00099
00100
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 }