00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooPolynomial.h"
00021 #include "RooAddPdf.h"
00022 #include "RooProdPdf.h"
00023 #include "TCanvas.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026
00027
00028 class TestBasic316 : public RooFitTestUnit
00029 {
00030 public:
00031 TestBasic316(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Likelihood ratio projection plot",refFile,writeRef,verbose) {} ;
00032 Bool_t testCode() {
00033
00034
00035
00036
00037
00038 RooRealVar x("x","x",-5,5) ;
00039 RooRealVar y("y","y",-5,5) ;
00040 RooRealVar z("z","z",-5,5) ;
00041
00042
00043 RooGaussian gx("gx","gx",x,RooConst(0),RooConst(1)) ;
00044 RooGaussian gy("gy","gy",y,RooConst(0),RooConst(1)) ;
00045 RooGaussian gz("gz","gz",z,RooConst(0),RooConst(1)) ;
00046 RooProdPdf sig("sig","sig",RooArgSet(gx,gy,gz)) ;
00047
00048
00049 RooPolynomial px("px","px",x,RooArgSet(RooConst(-0.1),RooConst(0.004))) ;
00050 RooPolynomial py("py","py",y,RooArgSet(RooConst(0.1),RooConst(-0.004))) ;
00051 RooPolynomial pz("pz","pz",z) ;
00052 RooProdPdf bkg("bkg","bkg",RooArgSet(px,py,pz)) ;
00053
00054
00055 RooRealVar fsig("fsig","signal fraction",0.1,0.,1.) ;
00056 RooAddPdf model("model","model",RooArgList(sig,bkg),fsig) ;
00057
00058 RooDataSet* data = model.generate(RooArgSet(x,y,z),20000) ;
00059
00060
00061
00062
00063
00064
00065
00066 RooPlot* frame = x.frame(Title("Projection of 3D data and pdf on X"),Bins(40)) ;
00067 data->plotOn(frame) ;
00068 model.plotOn(frame) ;
00069
00070
00071
00072
00073
00074
00075
00076
00077 RooAbsPdf* sigyz = sig.createProjection(x) ;
00078 RooAbsPdf* totyz = model.createProjection(x) ;
00079
00080
00081 RooFormulaVar llratio_func("llratio","log10(@0)-log10(@1)",RooArgList(*sigyz,*totyz)) ;
00082
00083
00084
00085
00086
00087
00088
00089 data->addColumn(llratio_func) ;
00090
00091
00092 RooDataSet* dataSel = (RooDataSet*) data->reduce(Cut("llratio>0.7")) ;
00093
00094
00095 RooPlot* frame2 = x.frame(Title("Same projection on X with LLratio(y,z)>0.7"),Bins(40)) ;
00096
00097
00098 dataSel->plotOn(frame2) ;
00099
00100
00101
00102
00103
00104
00105
00106 RooDataSet* mcprojData = model.generate(RooArgSet(x,y,z),10000) ;
00107
00108
00109 mcprojData->addColumn(llratio_func) ;
00110 RooDataSet* mcprojDataSel = (RooDataSet*) mcprojData->reduce(Cut("llratio>0.7")) ;
00111
00112
00113
00114 model.plotOn(frame2,ProjWData(*mcprojDataSel)) ;
00115
00116
00117 regPlot(frame,"rf316_plot1") ;
00118 regPlot(frame2,"rf316_plot2") ;
00119
00120 delete data ;
00121 delete dataSel ;
00122 delete mcprojData ;
00123 delete mcprojDataSel ;
00124 delete sigyz ;
00125 delete totyz ;
00126
00127 return kTRUE ;
00128 }
00129 } ;