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 "RooCategory.h"
00019 #include "RooEfficiency.h"
00020 #include "RooPolynomial.h"
00021 #include "RooProdPdf.h"
00022 #include "RooFormulaVar.h"
00023 #include "TCanvas.h"
00024 #include "TH1.h"
00025 #include "RooPlot.h"
00026 using namespace RooFit ;
00027
00028
00029
00030 class TestBasic702 : public RooFitTestUnit
00031 {
00032 public:
00033 TestBasic702(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Efficiency operator p.d.f. 2D",refFile,writeRef,verbose) {} ;
00034 Bool_t testCode() {
00035
00036 Bool_t flat=kFALSE ;
00037
00038
00039
00040
00041
00042 RooRealVar x("x","x",-10,10) ;
00043 RooRealVar y("y","y",-10,10) ;
00044
00045
00046 RooRealVar ax("ax","ay",0.6,0,1) ;
00047 RooRealVar bx("bx","by",5) ;
00048 RooRealVar cx("cx","cy",-1,-10,10) ;
00049
00050 RooRealVar ay("ay","ay",0.2,0,1) ;
00051 RooRealVar by("by","by",5) ;
00052 RooRealVar cy("cy","cy",-1,-10,10) ;
00053
00054 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)) ;
00055
00056
00057 RooCategory cut("cut","cutr") ;
00058 cut.defineType("accept",1) ;
00059 cut.defineType("reject",0) ;
00060
00061
00062
00063
00064
00065
00066
00067 RooEfficiency effPdf("effPdf","effPdf",effFunc,cut,"accept") ;
00068
00069
00070
00071
00072
00073
00074
00075
00076 RooPolynomial shapePdfX("shapePdfX","shapePdfX",x,RooConst(flat?0:-0.095)) ;
00077 RooPolynomial shapePdfY("shapePdfY","shapePdfY",y,RooConst(flat?0:+0.095)) ;
00078 RooProdPdf shapePdf("shapePdf","shapePdf",RooArgSet(shapePdfX,shapePdfY)) ;
00079 RooProdPdf model("model","model",shapePdf,Conditional(effPdf,cut)) ;
00080
00081
00082 RooDataSet* data = model.generate(RooArgSet(x,y,cut),10000) ;
00083
00084
00085
00086
00087
00088
00089
00090 effPdf.fitTo(*data,ConditionalObservables(RooArgSet(x,y))) ;
00091
00092
00093
00094
00095
00096
00097
00098 TH1* hh_data_all = data->createHistogram("hh_data_all",x,Binning(8),YVar(y,Binning(8))) ;
00099 TH1* hh_data_sel = data->createHistogram("hh_data_sel",x,Binning(8),YVar(y,Binning(8)),Cut("cut==cut::accept")) ;
00100 TH1* hh_eff = effFunc.createHistogram("hh_eff",x,Binning(50),YVar(y,Binning(50))) ;
00101
00102
00103 hh_data_all->SetMinimum(0) ;
00104 hh_data_sel->SetMinimum(0) ;
00105 hh_eff->SetMinimum(0) ;
00106 hh_eff->SetLineColor(kBlue) ;
00107
00108 regTH(hh_data_all,"rf702_hh_data_all") ;
00109 regTH(hh_data_sel,"rf702_hh_data_sel") ;
00110 regTH(hh_eff,"rf702_hh_eff") ;
00111
00112 delete data ;
00113
00114 return kTRUE;
00115
00116 }
00117 } ;