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 "RooExponential.h"
00020 #include "RooEffProd.h"
00021 #include "RooFormulaVar.h"
00022 #include "TCanvas.h"
00023 #include "RooPlot.h"
00024 using namespace RooFit ;
00025
00026
00027
00028 class TestBasic703 : public RooFitTestUnit
00029 {
00030 public:
00031 TestBasic703(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Efficiency product operator p.d.f",refFile,writeRef,verbose) {} ;
00032 Bool_t testCode() {
00033
00034
00035
00036
00037
00038 RooRealVar t("t","t",0,5) ;
00039
00040
00041 RooRealVar tau("tau","tau",-1.54,-4,-0.1) ;
00042 RooExponential model("model","model",t,tau) ;
00043
00044
00045
00046
00047
00048
00049
00050 RooFormulaVar eff("eff","0.5*(TMath::Erf((t-1)/0.5)+1)",t) ;
00051
00052
00053
00054
00055
00056
00057
00058 RooEffProd modelEff("modelEff","model with efficiency",model,eff) ;
00059
00060
00061
00062
00063
00064
00065 RooPlot* frame1 = t.frame(Title("Efficiency")) ;
00066 eff.plotOn(frame1,LineColor(kRed)) ;
00067
00068 RooPlot* frame2 = t.frame(Title("Pdf with and without efficiency")) ;
00069
00070 model.plotOn(frame2,LineStyle(kDashed)) ;
00071 modelEff.plotOn(frame2) ;
00072
00073
00074
00075
00076
00077
00078
00079
00080 RooDataSet* data = modelEff.generate(t,10000) ;
00081
00082
00083 modelEff.fitTo(*data) ;
00084
00085
00086 RooPlot* frame3 = t.frame(Title("Fitted pdf with efficiency")) ;
00087 data->plotOn(frame3) ;
00088 modelEff.plotOn(frame3) ;
00089
00090
00091 regPlot(frame1,"rf703_plot1") ;
00092 regPlot(frame2,"rf703_plot2") ;
00093 regPlot(frame3,"rf703_plot3") ;
00094
00095
00096 delete data ;
00097 return kTRUE ;
00098 }
00099 } ;