rf310_sliceplot.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #309
00004 // 
00005 // Projecting p.d.f and data slices in discrete observables 
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussModel.h"
00019 #include "RooDecay.h"
00020 #include "RooBMixDecay.h"
00021 #include "RooCategory.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026 
00027 
00028 void rf310_sliceplot()
00029 {
00030 
00031 
00032   // C r e a t e   B   d e c a y   p d f   w it h   m i x i n g 
00033   // ----------------------------------------------------------
00034 
00035   // Decay time observables
00036   RooRealVar dt("dt","dt",-20,20) ;
00037 
00038   // Discrete observables mixState (B0tag==B0reco?) and tagFlav (B0tag==B0(bar)?)
00039   RooCategory mixState("mixState","B0/B0bar mixing state") ;
00040   RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
00041 
00042   // Define state labels of discrete observables
00043   mixState.defineType("mixed",-1) ;
00044   mixState.defineType("unmixed",1) ;
00045   tagFlav.defineType("B0",1) ;
00046   tagFlav.defineType("B0bar",-1) ;
00047 
00048   // Model parameters
00049   RooRealVar dm("dm","delta m(B)",0.472,0.,1.0) ;
00050   RooRealVar tau("tau","B0 decay time",1.547,1.0,2.0) ;
00051   RooRealVar w("w","Flavor Mistag rate",0.03,0.0,1.0) ;
00052   RooRealVar dw("dw","Flavor Mistag rate difference between B0 and B0bar",0.01) ;
00053 
00054   // Build a gaussian resolution model
00055   RooRealVar bias1("bias1","bias1",0) ;
00056   RooRealVar sigma1("sigma1","sigma1",0.01) ;  
00057   RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
00058 
00059   // Construct a decay pdf, smeared with single gaussian resolution model
00060   RooBMixDecay bmix_gm1("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,RooBMixDecay::DoubleSided) ;
00061   
00062   // Generate BMixing data with above set of event errors
00063   RooDataSet *data = bmix_gm1.generate(RooArgSet(dt,tagFlav,mixState),20000) ;
00064 
00065 
00066 
00067   // P l o t   f u l l   d e c a y   d i s t r i b u t i o n 
00068   // ----------------------------------------------------------
00069 
00070   // Create frame, plot data and pdf projection (integrated over tagFlav and mixState)
00071   RooPlot* frame = dt.frame(Title("Inclusive decay distribution")) ;
00072   data->plotOn(frame) ;
00073   bmix_gm1.plotOn(frame) ;
00074 
00075 
00076 
00077   // P l o t   d e c a y   d i s t r .   f o r   m i x e d   a n d   u n m i x e d   s l i c e   o f   m i x S t a t e
00078   // ------------------------------------------------------------------------------------------------------------------
00079 
00080   // Create frame, plot data (mixed only)
00081   RooPlot* frame2 = dt.frame(Title("Decay distribution of mixed events")) ;
00082   data->plotOn(frame2,Cut("mixState==mixState::mixed")) ;
00083 
00084   // Position slice in mixState at "mixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
00085   bmix_gm1.plotOn(frame2,Slice(mixState,"mixed")) ;
00086 
00087   // Create frame, plot data (unmixed only)
00088   RooPlot* frame3 = dt.frame(Title("Decay distribution of unmixed events")) ;
00089   data->plotOn(frame3,Cut("mixState==mixState::unmixed")) ;
00090 
00091   // Position slice in mixState at "unmixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
00092   bmix_gm1.plotOn(frame3,Slice(mixState,"unmixed")) ;
00093 
00094 
00095 
00096   TCanvas* c = new TCanvas("rf310_sliceplot","rf310_sliceplot",1200,400) ;
00097   c->Divide(3) ;
00098   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame->Draw() ;
00099   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame2->Draw() ;
00100   c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame3->Draw() ;
00101 
00102 
00103 }

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