00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooGaussian.h"
00018 #include "RooConstVar.h"
00019 #include "RooCategory.h"
00020 #include "RooEfficiency.h"
00021 #include "RooPolynomial.h"
00022 #include "RooProdPdf.h"
00023 #include "RooFormulaVar.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 rf702_efficiencyfit_2D(Bool_t flat=kFALSE)
00032 {
00033
00034
00035
00036
00037 RooRealVar x("x","x",-10,10) ;
00038 RooRealVar y("y","y",-10,10) ;
00039
00040
00041 RooRealVar ax("ax","ay",0.6,0,1) ;
00042 RooRealVar bx("bx","by",5) ;
00043 RooRealVar cx("cx","cy",-1,-10,10) ;
00044
00045 RooRealVar ay("ay","ay",0.2,0,1) ;
00046 RooRealVar by("by","by",5) ;
00047 RooRealVar cy("cy","cy",-1,-10,10) ;
00048
00049 RooFormulaVar effFunc("effFunc","((1-ax)+ax*cos((x-cx)/bx))*((1-ay)+ay*cos((y-cy)/by))",RooArgList(ax,bx,cx,x,ay,by,cy,y)) ;
00050
00051
00052 RooCategory cut("cut","cutr") ;
00053 cut.defineType("accept",1) ;
00054 cut.defineType("reject",0) ;
00055
00056
00057
00058
00059
00060
00061
00062 RooEfficiency effPdf("effPdf","effPdf",effFunc,cut,"accept") ;
00063
00064
00065
00066
00067
00068
00069
00070
00071 RooPolynomial shapePdfX("shapePdfX","shapePdfX",x,RooConst(flat?0:-0.095)) ;
00072 RooPolynomial shapePdfY("shapePdfY","shapePdfY",y,RooConst(flat?0:+0.095)) ;
00073 RooProdPdf shapePdf("shapePdf","shapePdf",RooArgSet(shapePdfX,shapePdfY)) ;
00074 RooProdPdf model("model","model",shapePdf,Conditional(effPdf,cut)) ;
00075
00076
00077 RooDataSet* data = model.generate(RooArgSet(x,y,cut),10000) ;
00078
00079
00080
00081
00082
00083
00084
00085 effPdf.fitTo(*data,ConditionalObservables(RooArgSet(x,y))) ;
00086
00087
00088
00089
00090
00091
00092
00093 TH1* hh_data_all = data->createHistogram("hh_data_all",x,Binning(8),YVar(y,Binning(8))) ;
00094 TH1* hh_data_sel = data->createHistogram("hh_data_sel",x,Binning(8),YVar(y,Binning(8)),Cut("cut==cut::accept")) ;
00095 TH1* hh_eff = effFunc.createHistogram("hh_eff",x,Binning(50),YVar(y,Binning(50))) ;
00096
00097
00098 hh_data_all->SetMinimum(0) ;
00099 hh_data_sel->SetMinimum(0) ;
00100 hh_eff->SetMinimum(0) ;
00101 hh_eff->SetLineColor(kBlue) ;
00102
00103
00104
00105
00106 TCanvas* ca = new TCanvas("rf702_efficiency_2D","rf702_efficiency_2D",1200,400) ;
00107 ca->Divide(3) ;
00108 ca->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_data_all->GetZaxis()->SetTitleOffset(1.8) ; hh_data_all->Draw("lego") ;
00109 ca->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_data_sel->GetZaxis()->SetTitleOffset(1.8) ; hh_data_sel->Draw("lego") ;
00110 ca->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_eff->GetZaxis()->SetTitleOffset(1.8) ; hh_eff->Draw("surf") ;
00111
00112 return ;
00113
00114
00115 }