ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf708_bphysics.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'SPECIAL PDFS' RooFit tutorial macro #708
4 ///
5 /// Special decay pdf for B physics with mixing and/or CP violation
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 "RooGaussian.h"
19 #include "RooConstVar.h"
20 #include "RooCategory.h"
21 #include "RooBMixDecay.h"
22 #include "RooBCPEffDecay.h"
23 #include "RooBDecay.h"
24 #include "RooFormulaVar.h"
25 #include "RooTruthModel.h"
26 #include "TCanvas.h"
27 #include "TAxis.h"
28 #include "RooPlot.h"
29 using namespace RooFit ;
30 
31 void rf708_bphysics()
32 {
33  ////////////////////////////////////////////////////
34  // B - D e c a y w i t h m i x i n g //
35  ////////////////////////////////////////////////////
36 
37  // C o n s t r u c t p d f
38  // -------------------------
39 
40  // Observable
41  RooRealVar dt("dt","dt",-10,10) ;
42  dt.setBins(40) ;
43 
44  // Parameters
45  RooRealVar dm("dm","delta m(B0)",0.472) ;
46  RooRealVar tau("tau","tau (B0)",1.547) ;
47  RooRealVar w("w","flavour mistag rate",0.1) ;
48  RooRealVar dw("dw","delta mistag rate for B0/B0bar",0.1) ;
49 
50  RooCategory mixState("mixState","B0/B0bar mixing state") ;
51  mixState.defineType("mixed",-1) ;
52  mixState.defineType("unmixed",1) ;
53 
54  RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
55  tagFlav.defineType("B0",1) ;
56  tagFlav.defineType("B0bar",-1) ;
57 
58  // Use delta function resolution model
59  RooTruthModel tm("tm","truth model",dt) ;
60 
61  // Construct Bdecay with mixing
62  RooBMixDecay bmix("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,tm,RooBMixDecay::DoubleSided) ;
63 
64 
65 
66  // P l o t p d f i n v a r i o u s s l i c e s
67  // ---------------------------------------------------
68 
69  // Generate some data
70  RooDataSet* data = bmix.generate(RooArgSet(dt,mixState,tagFlav),10000) ;
71 
72  // Plot B0 and B0bar tagged data separately
73  // For all plots below B0 and B0 tagged data will look somewhat differently
74  // if the flavor tagging mistag rate for B0 and B0 is different (i.e. dw!=0)
75  RooPlot* frame1 = dt.frame(Title("B decay distribution with mixing (B0/B0bar)")) ;
76 
77  data->plotOn(frame1,Cut("tagFlav==tagFlav::B0")) ;
78  bmix.plotOn(frame1,Slice(tagFlav,"B0")) ;
79 
80  data->plotOn(frame1,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
81  bmix.plotOn(frame1,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
82 
83 
84  // Plot mixed slice for B0 and B0bar tagged data separately
85  RooPlot* frame2 = dt.frame(Title("B decay distribution of mixed events (B0/B0bar)")) ;
86 
87  data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0")) ;
88  bmix.plotOn(frame2,Slice(tagFlav,"B0"),Slice(mixState,"mixed")) ;
89 
90  data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
91  bmix.plotOn(frame2,Slice(tagFlav,"B0bar"),Slice(mixState,"mixed"),LineColor(kCyan)) ;
92 
93 
94  // Plot unmixed slice for B0 and B0bar tagged data separately
95  RooPlot* frame3 = dt.frame(Title("B decay distribution of unmixed events (B0/B0bar)")) ;
96 
97  data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0")) ;
98  bmix.plotOn(frame3,Slice(tagFlav,"B0"),Slice(mixState,"unmixed")) ;
99 
100  data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
101  bmix.plotOn(frame3,Slice(tagFlav,"B0bar"),Slice(mixState,"unmixed"),LineColor(kCyan)) ;
102 
103 
104 
105 
106 
107  ///////////////////////////////////////////////////////
108  // B - D e c a y w i t h C P v i o l a t i o n //
109  ///////////////////////////////////////////////////////
110 
111  // C o n s t r u c t p d f
112  // -------------------------
113 
114  // Additional parameters needed for B decay with CPV
115  RooRealVar CPeigen("CPeigen","CP eigen value",-1) ;
116  RooRealVar absLambda("absLambda","|lambda|",1,0,2) ;
117  RooRealVar argLambda("absLambda","|lambda|",0.7,-1,1) ;
118  RooRealVar effR("effR","B0/B0bar reco efficiency ratio",1) ;
119 
120  // Construct Bdecay with CP violation
121  RooBCPEffDecay bcp("bcp","bcp", dt, tagFlav, tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, tm, RooBCPEffDecay::DoubleSided) ;
122 
123 
124 
125  // P l o t s c e n a r i o 1 - s i n ( 2 b ) = 0 . 7 , | l | = 1
126  // ---------------------------------------------------------------------------
127 
128  // Generate some data
129  RooDataSet* data2 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
130 
131  // Plot B0 and B0bar tagged data separately
132  RooPlot* frame4 = dt.frame(Title("B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)")) ;
133 
134  data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0")) ;
135  bcp.plotOn(frame4,Slice(tagFlav,"B0")) ;
136 
137  data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
138  bcp.plotOn(frame4,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
139 
140 
141 
142  // P l o t s c e n a r i o 2 - s i n ( 2 b ) = 0 . 7 , | l | = 0 . 7
143  // -------------------------------------------------------------------------------
144 
145  absLambda=0.7 ;
146 
147  // Generate some data
148  RooDataSet* data3 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
149 
150  // Plot B0 and B0bar tagged data separately (sin2b = 0.7 plus direct CPV |l|=0.5)
151  RooPlot* frame5 = dt.frame(Title("B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)")) ;
152 
153  data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0")) ;
154  bcp.plotOn(frame5,Slice(tagFlav,"B0")) ;
155 
156  data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
157  bcp.plotOn(frame5,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
158 
159 
160 
161  //////////////////////////////////////////////////////////////////////////////////
162  // G e n e r i c B d e c a y w i t h u s e r c o e f f i c i e n t s //
163  //////////////////////////////////////////////////////////////////////////////////
164 
165  // C o n s t r u c t p d f
166  // -------------------------
167 
168  // Model parameters
169  RooRealVar DGbG("DGbG","DGamma/GammaAvg",0.5,-1,1);
170  RooRealVar Adir("Adir","-[1-abs(l)**2]/[1+abs(l)**2]",0);
171  RooRealVar Amix("Amix","2Im(l)/[1+abs(l)**2]",0.7);
172  RooRealVar Adel("Adel","2Re(l)/[1+abs(l)**2]",0.7);
173 
174  // Derived input parameters for pdf
175  RooFormulaVar DG("DG","Delta Gamma","@1/@0",RooArgList(tau,DGbG));
176 
177  // Construct coefficient functions for sin,cos,sinh modulations of decay distribution
178  RooFormulaVar fsin("fsin","fsin","@0*@1*(1-2*@2)",RooArgList(Amix,tagFlav,w));
179  RooFormulaVar fcos("fcos","fcos","@0*@1*(1-2*@2)",RooArgList(Adir,tagFlav,w));
180  RooFormulaVar fsinh("fsinh","fsinh","@0",RooArgList(Adel));
181 
182  // Construct generic B decay pdf using above user coefficients
183  RooBDecay bcpg("bcpg","bcpg",dt,tau,DG,RooConst(1),fsinh,fcos,fsin,dm,tm,RooBDecay::DoubleSided);
184 
185 
186 
187  // P l o t - I m ( l ) = 0 . 7 , R e ( l ) = 0 . 7 | l | = 1, d G / G = 0 . 5
188  // -------------------------------------------------------------------------------------
189 
190  // Generate some data
191  RooDataSet* data4 = bcpg.generate(RooArgSet(dt,tagFlav),10000) ;
192 
193  // Plot B0 and B0bar tagged data separately
194  RooPlot* frame6 = dt.frame(Title("B decay distribution with CPV(Im(l)=0.7,Re(l)=0.7,|l|=1,dG/G=0.5) (B0/B0bar)")) ;
195 
196  data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0")) ;
197  bcpg.plotOn(frame6,Slice(tagFlav,"B0")) ;
198 
199  data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
200  bcpg.plotOn(frame6,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
201 
202 
203 
204 
205  TCanvas* c = new TCanvas("rf708_bphysics","rf708_bphysics",1200,800) ;
206  c->Divide(3,2) ;
207  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.6) ; frame1->Draw() ;
208  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
209  c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
210  c->cd(4) ; gPad->SetLeftMargin(0.15) ; frame4->GetYaxis()->SetTitleOffset(1.6) ; frame4->Draw() ;
211  c->cd(5) ; gPad->SetLeftMargin(0.15) ; frame5->GetYaxis()->SetTitleOffset(1.6) ; frame5->Draw() ;
212  c->cd(6) ; gPad->SetLeftMargin(0.15) ; frame6->GetYaxis()->SetTitleOffset(1.6) ; frame6->Draw() ;
213 
214 }
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
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineColor(Color_t color)
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
PDF describing decay time distribution of B meson including effects of standard model CP violation...
RooCmdArg Title(const char *name)
Most general description of B decay time distribution with effects of CP violation, mixing and life time differences.
Definition: RooBDecay.h:24
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
Definition: Rtypes.h:61
RooConstVar & RooConst(Double_t val)
#define gPad
Definition: TVirtualPad.h:288
RooCmdArg MarkerColor(Color_t color)
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