00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 #ifndef __CINT__
00019 #include "RooGlobalFunc.h"
00020 #endif
00021 #include "RooRealVar.h"
00022 #include "RooDataSet.h"
00023 #include "RooGaussModel.h"
00024 #include "RooAddModel.h"
00025 #include "RooTruthModel.h"
00026 #include "RooDecay.h"
00027 #include "RooPlot.h"
00028 #include "TCanvas.h"
00029 #include "TAxis.h"
00030 #include "TH1.h"
00031 using namespace RooFit ;
00032
00033
00034
00035 void rf209_anaconv()
00036 {
00037
00038
00039
00040
00041 RooRealVar dt("dt","dt",-10,10) ;
00042 RooRealVar tau("tau","tau",1.548) ;
00043
00044
00045 RooTruthModel tm("tm","truth model",dt) ;
00046
00047
00048 RooDecay decay_tm("decay_tm","decay",dt,tau,tm,RooDecay::DoubleSided) ;
00049
00050
00051 RooPlot* frame = dt.frame(Title("Bdecay (x) resolution")) ;
00052 decay_tm.plotOn(frame,LineStyle(kDashed)) ;
00053
00054
00055
00056
00057
00058
00059 RooRealVar bias1("bias1","bias1",0) ;
00060 RooRealVar sigma1("sigma1","sigma1",1) ;
00061 RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
00062
00063
00064 RooDecay decay_gm1("decay_gm1","decay",dt,tau,gm1,RooDecay::DoubleSided) ;
00065
00066
00067 decay_gm1.plotOn(frame) ;
00068
00069
00070
00071
00072
00073
00074 RooRealVar bias2("bias2","bias2",0) ;
00075 RooRealVar sigma2("sigma2","sigma2",5) ;
00076 RooGaussModel gm2("gm2","gauss model 2",dt,bias2,sigma2) ;
00077
00078
00079 RooRealVar gm1frac("gm1frac","fraction of gm1",0.5) ;
00080 RooAddModel gmsum("gmsum","sum of gm1 and gm2",RooArgList(gm1,gm2),gm1frac) ;
00081
00082
00083 RooDecay decay_gmsum("decay_gmsum","decay",dt,tau,gmsum,RooDecay::DoubleSided) ;
00084
00085
00086 decay_gmsum.plotOn(frame,LineColor(kRed)) ;
00087
00088
00089
00090 new TCanvas("rf209_anaconv","rf209_anaconv",600, 600);
00091 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00092
00093 }