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 "RooGaussModel.h"
00021 #include "RooDecay.h"
00022 #include "RooLandau.h"
00023 #include "RooPlot.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "TH2D.h"
00027 using namespace RooFit ;
00028
00029
00030
00031 void rf306_condpereventerrors()
00032 {
00033
00034
00035
00036
00037 RooRealVar dt("dt","dt",-10,10) ;
00038 RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
00039
00040
00041 RooRealVar bias("bias","bias",0,-10,10) ;
00042 RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
00043 RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
00044
00045
00046 RooRealVar tau("tau","tau",1.548) ;
00047 RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
00048
00049
00050
00051
00052
00053
00054
00055 RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
00056 RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
00057
00058
00059
00060
00061
00062
00063
00064 RooDataSet* data = decay_gm.generate(dt,ProtoData(*expDataDterr)) ;
00065
00066
00067
00068
00069
00070
00071
00072 decay_gm.fitTo(*data,ConditionalObservables(dterr)) ;
00073
00074
00075
00076
00077
00078
00079
00080
00081 TH1* hh_decay = decay_gm.createHistogram("hh_decay",dt,Binning(50),YVar(dterr,Binning(50))) ;
00082 hh_decay->SetLineColor(kBlue) ;
00083
00084
00085
00086 RooPlot* frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr")) ;
00087 for (Int_t ibin=0 ; ibin<100 ; ibin+=20) {
00088 dterr.setBin(ibin) ;
00089 decay_gm.plotOn(frame,Normalization(5.)) ;
00090 }
00091
00092
00093
00094 RooPlot* frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt")) ;
00095 data->plotOn(frame2) ;
00096
00097
00098
00099
00100
00101
00102 decay_gm.plotOn(frame2,ProjWData(*expDataDterr,kTRUE)) ;
00103
00104
00105
00106
00107 TCanvas* c = new TCanvas("rf306_condpereventerrors","rf306_condperventerrors",1200, 400);
00108 c->Divide(3) ;
00109 c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_decay->GetZaxis()->SetTitleOffset(2.5) ; hh_decay->Draw("surf") ;
00110 c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00111 c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
00112
00113
00114 }
00115