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 "RooExponential.h"
00021 #include "RooEffProd.h"
00022 #include "RooFormulaVar.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 using namespace RooFit ;
00027
00028
00029 void rf703_effpdfprod()
00030 {
00031
00032
00033
00034
00035 RooRealVar t("t","t",0,5) ;
00036
00037
00038 RooRealVar tau("tau","tau",-1.54,-4,-0.1) ;
00039 RooExponential model("model","model",t,tau) ;
00040
00041
00042
00043
00044
00045
00046
00047 RooFormulaVar eff("eff","0.5*(TMath::Erf((t-1)/0.5)+1)",t) ;
00048
00049
00050
00051
00052
00053
00054
00055 RooEffProd modelEff("modelEff","model with efficiency",model,eff) ;
00056
00057
00058
00059
00060
00061
00062 RooPlot* frame1 = t.frame(Title("Efficiency")) ;
00063 eff.plotOn(frame1,LineColor(kRed)) ;
00064
00065 RooPlot* frame2 = t.frame(Title("Pdf with and without efficiency")) ;
00066
00067 model.plotOn(frame2,LineStyle(kDashed)) ;
00068 modelEff.plotOn(frame2) ;
00069
00070
00071
00072
00073
00074
00075
00076
00077 RooDataSet* data = modelEff.generate(t,10000) ;
00078
00079
00080 modelEff.fitTo(*data) ;
00081
00082
00083 RooPlot* frame3 = t.frame(Title("Fitted pdf with efficiency")) ;
00084 data->plotOn(frame3) ;
00085 modelEff.plotOn(frame3) ;
00086
00087
00088 TCanvas* c = new TCanvas("rf703_effpdfprod","rf703_effpdfprod",1200,400) ;
00089 c->Divide(3) ;
00090 c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
00091 c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
00092 c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
00093
00094 }