ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf310_sliceplot.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #309
4 ///
5 /// Projecting p.d.f and data slices in discrete observables
6 ///
7 ///
8 ///
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussModel.h"
19 #include "RooDecay.h"
20 #include "RooBMixDecay.h"
21 #include "RooCategory.h"
22 #include "TCanvas.h"
23 #include "TAxis.h"
24 #include "RooPlot.h"
25 using namespace RooFit ;
26 
27 
28 void rf310_sliceplot()
29 {
30 
31 
32  // C r e a t e B d e c a y p d f w it h m i x i n g
33  // ----------------------------------------------------------
34 
35  // Decay time observables
36  RooRealVar dt("dt","dt",-20,20) ;
37 
38  // Discrete observables mixState (B0tag==B0reco?) and tagFlav (B0tag==B0(bar)?)
39  RooCategory mixState("mixState","B0/B0bar mixing state") ;
40  RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
41 
42  // Define state labels of discrete observables
43  mixState.defineType("mixed",-1) ;
44  mixState.defineType("unmixed",1) ;
45  tagFlav.defineType("B0",1) ;
46  tagFlav.defineType("B0bar",-1) ;
47 
48  // Model parameters
49  RooRealVar dm("dm","delta m(B)",0.472,0.,1.0) ;
50  RooRealVar tau("tau","B0 decay time",1.547,1.0,2.0) ;
51  RooRealVar w("w","Flavor Mistag rate",0.03,0.0,1.0) ;
52  RooRealVar dw("dw","Flavor Mistag rate difference between B0 and B0bar",0.01) ;
53 
54  // Build a gaussian resolution model
55  RooRealVar bias1("bias1","bias1",0) ;
56  RooRealVar sigma1("sigma1","sigma1",0.01) ;
57  RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
58 
59  // Construct a decay pdf, smeared with single gaussian resolution model
60  RooBMixDecay bmix_gm1("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,RooBMixDecay::DoubleSided) ;
61 
62  // Generate BMixing data with above set of event errors
63  RooDataSet *data = bmix_gm1.generate(RooArgSet(dt,tagFlav,mixState),20000) ;
64 
65 
66 
67  // 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
68  // ----------------------------------------------------------
69 
70  // Create frame, plot data and pdf projection (integrated over tagFlav and mixState)
71  RooPlot* frame = dt.frame(Title("Inclusive decay distribution")) ;
72  data->plotOn(frame) ;
73  bmix_gm1.plotOn(frame) ;
74 
75 
76 
77  // 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
78  // ------------------------------------------------------------------------------------------------------------------
79 
80  // Create frame, plot data (mixed only)
81  RooPlot* frame2 = dt.frame(Title("Decay distribution of mixed events")) ;
82  data->plotOn(frame2,Cut("mixState==mixState::mixed")) ;
83 
84  // Position slice in mixState at "mixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
85  bmix_gm1.plotOn(frame2,Slice(mixState,"mixed")) ;
86 
87  // Create frame, plot data (unmixed only)
88  RooPlot* frame3 = dt.frame(Title("Decay distribution of unmixed events")) ;
89  data->plotOn(frame3,Cut("mixState==mixState::unmixed")) ;
90 
91  // Position slice in mixState at "unmixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
92  bmix_gm1.plotOn(frame3,Slice(mixState,"unmixed")) ;
93 
94 
95 
96  TCanvas* c = new TCanvas("rf310_sliceplot","rf310_sliceplot",1200,400) ;
97  c->Divide(3) ;
98  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame->Draw() ;
99  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame2->Draw() ;
100  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.4) ; gPad->SetLogy() ; frame3->Draw() ;
101 
102 
103 }
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:245
RooCmdArg Cut(const char *cutSpec)
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
RooCmdArg Title(const char *name)
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Definition: RooBMixDecay.h:23
friend class RooArgSet
Definition: RooAbsArg.h:469
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
Definition: RooAbsData.cxx:626
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
tuple w
Definition: qtexample.py:51
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
#define gPad
Definition: TVirtualPad.h:288
RooCmdArg Slice(const RooArgSet &sliceSet)
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26