Logo ROOT   6.10/09
Reference Guide
rf310_sliceplot.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #309
4 //
5 // Projecting p.d.f and data slices in discrete observables
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
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 "RooPlot.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic310 : public RooFitTestUnit
28 {
29 public:
30  TestBasic310(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Data and p.d.f projection in category slice",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // C r e a t e B d e c a y p d f w it h m i x i n g
34  // ----------------------------------------------------------
35 
36  // Decay time observables
37  RooRealVar dt("dt","dt",-20,20) ;
38 
39  // Discrete observables mixState (B0tag==B0reco?) and tagFlav (B0tag==B0(bar)?)
40  RooCategory mixState("mixState","B0/B0bar mixing state") ;
41  RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
42 
43  // Define state labels of discrete observables
44  mixState.defineType("mixed",-1) ;
45  mixState.defineType("unmixed",1) ;
46  tagFlav.defineType("B0",1) ;
47  tagFlav.defineType("B0bar",-1) ;
48 
49  // Model parameters
50  RooRealVar dm("dm","delta m(B)",0.472,0.,1.0) ;
51  RooRealVar tau("tau","B0 decay time",1.547,1.0,2.0) ;
52  RooRealVar w("w","Flavor Mistag rate",0.03,0.0,1.0) ;
53  RooRealVar dw("dw","Flavor Mistag rate difference between B0 and B0bar",0.01) ;
54 
55  // Build a gaussian resolution model
56  RooRealVar bias1("bias1","bias1",0) ;
57  RooRealVar sigma1("sigma1","sigma1",0.01) ;
58  RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
59 
60  // Construct a decay pdf, smeared with single gaussian resolution model
61  RooBMixDecay bmix_gm1("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,gm1,RooBMixDecay::DoubleSided) ;
62 
63  // Generate BMixing data with above set of event errors
64  RooDataSet *data = bmix_gm1.generate(RooArgSet(dt,tagFlav,mixState),20000) ;
65 
66 
67 
68  // 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
69  // ----------------------------------------------------------
70 
71  // Create frame, plot data and pdf projection (integrated over tagFlav and mixState)
72  RooPlot* frame = dt.frame(Title("Inclusive decay distribution")) ;
73  data->plotOn(frame) ;
74  bmix_gm1.plotOn(frame) ;
75 
76 
77 
78  // 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
79  // ------------------------------------------------------------------------------------------------------------------
80 
81  // Create frame, plot data (mixed only)
82  RooPlot* frame2 = dt.frame(Title("Decay distribution of mixed events")) ;
83  data->plotOn(frame2,Cut("mixState==mixState::mixed")) ;
84 
85  // Position slice in mixState at "mixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
86  bmix_gm1.plotOn(frame2,Slice(mixState,"mixed")) ;
87 
88  // Create frame, plot data (unmixed only)
89  RooPlot* frame3 = dt.frame(Title("Decay distribution of unmixed events")) ;
90  data->plotOn(frame3,Cut("mixState==mixState::unmixed")) ;
91 
92  // Position slice in mixState at "unmixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
93  bmix_gm1.plotOn(frame3,Slice(mixState,"unmixed")) ;
94 
95 
96  regPlot(frame,"rf310_plot1") ;
97  regPlot(frame2,"rf310_plot2") ;
98  regPlot(frame3,"rf310_plot3") ;
99 
100  delete data ;
101 
102  return kTRUE ;
103  }
104 } ;
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:552
RooCmdArg Cut(const char *cutSpec)
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Definition: RooBMixDecay.h:23
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
bool verbose
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:24
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg Slice(const RooArgSet &sliceSet)
const Bool_t kTRUE
Definition: RtypesCore.h:91
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26