rf209_anaconv.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #209
00004 // 
00005 // Decay function p.d.fs with optional B physics 
00006 // effects (mixing and CP violation) that can be
00007 // analytically convolved with e.g. Gaussian resolution 
00008 // functions
00009 // 
00010 // pdf1 = decay(t,tau) (x) delta(t)
00011 // pdf2 = decay(t,tau) (x) gauss(t,m,s)
00012 // pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
00013 // 
00014 // 07/2008 - Wouter Verkerke 
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   // B - p h y s i c s   p d f   w i t h   t r u t h   r e s o l u t i o n
00038   // ---------------------------------------------------------------------
00039 
00040   // Variables of decay p.d.f.
00041   RooRealVar dt("dt","dt",-10,10) ;
00042   RooRealVar tau("tau","tau",1.548) ;
00043 
00044   // Build a truth resolution model (delta function)
00045   RooTruthModel tm("tm","truth model",dt) ;
00046 
00047   // Construct decay(t) (x) delta(t)
00048   RooDecay decay_tm("decay_tm","decay",dt,tau,tm,RooDecay::DoubleSided) ;
00049 
00050   // Plot p.d.f. (dashed)
00051   RooPlot* frame = dt.frame(Title("Bdecay (x) resolution")) ;
00052   decay_tm.plotOn(frame,LineStyle(kDashed)) ;
00053 
00054 
00055   // B - p h y s i c s   p d f   w i t h   G a u s s i a n   r e s o l u t i o n
00056   // ----------------------------------------------------------------------------
00057 
00058   // Build a gaussian resolution model
00059   RooRealVar bias1("bias1","bias1",0) ;
00060   RooRealVar sigma1("sigma1","sigma1",1) ;
00061   RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
00062 
00063   // Construct decay(t) (x) gauss1(t)
00064   RooDecay decay_gm1("decay_gm1","decay",dt,tau,gm1,RooDecay::DoubleSided) ;
00065 
00066   // Plot p.d.f. 
00067   decay_gm1.plotOn(frame) ;
00068 
00069 
00070   // B - p h y s i c s   p d f   w i t h   d o u b l e   G a u s s i a n   r e s o l u t i o n
00071   // ------------------------------------------------------------------------------------------
00072 
00073   // Build another gaussian resolution model
00074   RooRealVar bias2("bias2","bias2",0) ;
00075   RooRealVar sigma2("sigma2","sigma2",5) ;
00076   RooGaussModel gm2("gm2","gauss model 2",dt,bias2,sigma2) ;
00077 
00078   // Build a composite resolution model f*gm1+(1-f)*gm2
00079   RooRealVar gm1frac("gm1frac","fraction of gm1",0.5) ;
00080   RooAddModel gmsum("gmsum","sum of gm1 and gm2",RooArgList(gm1,gm2),gm1frac) ;
00081 
00082   // Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
00083   RooDecay decay_gmsum("decay_gmsum","decay",dt,tau,gmsum,RooDecay::DoubleSided) ;
00084 
00085   // Plot p.d.f. (red)
00086   decay_gmsum.plotOn(frame,LineColor(kRed)) ;
00087 
00088 
00089   // Draw all frames on canvas
00090   new TCanvas("rf209_anaconv","rf209_anaconv",600, 600);
00091   gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00092 
00093 }

Generated on Tue Jul 5 15:45:03 2011 for ROOT_528-00b_version by  doxygen 1.5.1