ROOT  6.06/09
Reference Guide
rf708_bphysics.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
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 // 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 "RooGaussian.h"
19 #include "RooCategory.h"
20 #include "RooBMixDecay.h"
21 #include "RooBCPEffDecay.h"
22 #include "RooBDecay.h"
23 #include "RooFormulaVar.h"
24 #include "RooTruthModel.h"
25 #include "TCanvas.h"
26 #include "RooPlot.h"
27 using namespace RooFit ;
28 
29 // Elementary operations on a gaussian PDF
30 class TestBasic708 : public RooFitTestUnit
31 {
32 public:
33  TestBasic708(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("B Physics p.d.f.s",refFile,writeRef,verbose) {} ;
34  Bool_t testCode() {
35 
36  ////////////////////////////////////////////////////
37  // B - D e c a y w i t h m i x i n g //
38  ////////////////////////////////////////////////////
39 
40  // C o n s t r u c t p d f
41  // -------------------------
42 
43  // Observable
44  RooRealVar dt("dt","dt",-10,10) ;
45  dt.setBins(40) ;
46 
47  // Parameters
48  RooRealVar dm("dm","delta m(B0)",0.472) ;
49  RooRealVar tau("tau","tau (B0)",1.547) ;
50  RooRealVar w("w","flavour mistag rate",0.1) ;
51  RooRealVar dw("dw","delta mistag rate for B0/B0bar",0.1) ;
52 
53  RooCategory mixState("mixState","B0/B0bar mixing state") ;
54  mixState.defineType("mixed",-1) ;
55  mixState.defineType("unmixed",1) ;
56 
57  RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
58  tagFlav.defineType("B0",1) ;
59  tagFlav.defineType("B0bar",-1) ;
60 
61  // Use delta function resolution model
62  RooTruthModel tm("tm","truth model",dt) ;
63 
64  // Construct Bdecay with mixing
65  RooBMixDecay bmix("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,tm,RooBMixDecay::DoubleSided) ;
66 
67 
68 
69  // P l o t p d f i n v a r i o u s s l i c e s
70  // ---------------------------------------------------
71 
72  // Generate some data
73  RooDataSet* data = bmix.generate(RooArgSet(dt,mixState,tagFlav),10000) ;
74 
75  // Plot B0 and B0bar tagged data separately
76  // For all plots below B0 and B0 tagged data will look somewhat differently
77  // if the flavor tagging mistag rate for B0 and B0 is different (i.e. dw!=0)
78  RooPlot* frame1 = dt.frame(Title("B decay distribution with mixing (B0/B0bar)")) ;
79 
80  data->plotOn(frame1,Cut("tagFlav==tagFlav::B0")) ;
81  bmix.plotOn(frame1,Slice(tagFlav,"B0")) ;
82 
83  data->plotOn(frame1,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
84  bmix.plotOn(frame1,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
85 
86 
87  // Plot mixed slice for B0 and B0bar tagged data separately
88  RooPlot* frame2 = dt.frame(Title("B decay distribution of mixed events (B0/B0bar)")) ;
89 
90  data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0")) ;
91  bmix.plotOn(frame2,Slice(tagFlav,"B0"),Slice(mixState,"mixed")) ;
92 
93  data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
94  bmix.plotOn(frame2,Slice(tagFlav,"B0bar"),Slice(mixState,"mixed"),LineColor(kCyan),Name("alt")) ;
95 
96 
97  // Plot unmixed slice for B0 and B0bar tagged data separately
98  RooPlot* frame3 = dt.frame(Title("B decay distribution of unmixed events (B0/B0bar)")) ;
99 
100  data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0")) ;
101  bmix.plotOn(frame3,Slice(tagFlav,"B0"),Slice(mixState,"unmixed")) ;
102 
103  data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
104  bmix.plotOn(frame3,Slice(tagFlav,"B0bar"),Slice(mixState,"unmixed"),LineColor(kCyan),Name("alt")) ;
105 
106 
107 
108 
109 
110  ///////////////////////////////////////////////////////
111  // B - D e c a y w i t h C P v i o l a t i o n //
112  ///////////////////////////////////////////////////////
113 
114  // C o n s t r u c t p d f
115  // -------------------------
116 
117  // Additional parameters needed for B decay with CPV
118  RooRealVar CPeigen("CPeigen","CP eigen value",-1) ;
119  RooRealVar absLambda("absLambda","|lambda|",1,0,2) ;
120  RooRealVar argLambda("absLambda","|lambda|",0.7,-1,1) ;
121  RooRealVar effR("effR","B0/B0bar reco efficiency ratio",1) ;
122 
123  // Construct Bdecay with CP violation
124  RooBCPEffDecay bcp("bcp","bcp", dt, tagFlav, tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, tm, RooBCPEffDecay::DoubleSided) ;
125 
126 
127 
128  // P l o t s c e n a r i o 1 - s i n ( 2 b ) = 0 . 7 , | l | = 1
129  // ---------------------------------------------------------------------------
130 
131  // Generate some data
132  RooDataSet* data2 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
133 
134  // Plot B0 and B0bar tagged data separately
135  RooPlot* frame4 = dt.frame(Title("B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)")) ;
136 
137  data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0")) ;
138  bcp.plotOn(frame4,Slice(tagFlav,"B0")) ;
139 
140  data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
141  bcp.plotOn(frame4,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
142 
143 
144 
145  // P l o t s c e n a r i o 2 - s i n ( 2 b ) = 0 . 7 , | l | = 0 . 7
146  // -------------------------------------------------------------------------------
147 
148  absLambda=0.7 ;
149 
150  // Generate some data
151  RooDataSet* data3 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
152 
153  // Plot B0 and B0bar tagged data separately (sin2b = 0.7 plus direct CPV |l|=0.5)
154  RooPlot* frame5 = dt.frame(Title("B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)")) ;
155 
156  data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0")) ;
157  bcp.plotOn(frame5,Slice(tagFlav,"B0")) ;
158 
159  data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
160  bcp.plotOn(frame5,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
161 
162 
163 
164  //////////////////////////////////////////////////////////////////////////////////
165  // 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 //
166  //////////////////////////////////////////////////////////////////////////////////
167 
168  // C o n s t r u c t p d f
169  // -------------------------
170 
171  // Model parameters
172  RooRealVar DGbG("DGbG","DGamma/GammaAvg",0.5,-1,1);
173  RooRealVar Adir("Adir","-[1-abs(l)**2]/[1+abs(l)**2]",0);
174  RooRealVar Amix("Amix","2Im(l)/[1+abs(l)**2]",0.7);
175  RooRealVar Adel("Adel","2Re(l)/[1+abs(l)**2]",0.7);
176 
177  // Derived input parameters for pdf
178  RooFormulaVar DG("DG","Delta Gamma","@1/@0",RooArgList(tau,DGbG));
179 
180  // Construct coefficient functions for sin,cos,sinh modulations of decay distribution
181  RooFormulaVar fsin("fsin","fsin","@0*@1*(1-2*@2)",RooArgList(Amix,tagFlav,w));
182  RooFormulaVar fcos("fcos","fcos","@0*@1*(1-2*@2)",RooArgList(Adir,tagFlav,w));
183  RooFormulaVar fsinh("fsinh","fsinh","@0",RooArgList(Adel));
184 
185  // Construct generic B decay pdf using above user coefficients
186  RooBDecay bcpg("bcpg","bcpg",dt,tau,DG,RooConst(1),fsinh,fcos,fsin,dm,tm,RooBDecay::DoubleSided);
187 
188 
189 
190  // P l o t - I m ( l ) = 0 . 7 , R e ( l ) = 0 . 7 | l | = 1, d G / G = 0 . 5
191  // -------------------------------------------------------------------------------------
192 
193  // Generate some data
194  RooDataSet* data4 = bcpg.generate(RooArgSet(dt,tagFlav),10000) ;
195 
196  // Plot B0 and B0bar tagged data separately
197  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)")) ;
198 
199  data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0")) ;
200  bcpg.plotOn(frame6,Slice(tagFlav,"B0")) ;
201 
202  data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
203  bcpg.plotOn(frame6,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
204 
205 
206  regPlot(frame1,"rf708_plot1") ;
207  regPlot(frame2,"rf708_plot2") ;
208  regPlot(frame3,"rf708_plot3") ;
209  regPlot(frame4,"rf708_plot4") ;
210  regPlot(frame5,"rf708_plot5") ;
211  regPlot(frame6,"rf708_plot6") ;
212 
213  delete data ;
214  delete data2 ;
215  delete data3 ;
216  delete data4 ;
217 
218  return kTRUE ;
219  }
220 } ;
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineColor(Color_t color)
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
PDF describing decay time distribution of B meson including effects of standard model CP violation...
RooCmdArg Title(const char *name)
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:624
RooCmdArg Name(const char *name)
bool verbose
Definition: Rtypes.h:61
RooConstVar & RooConst(Double_t val)
RooCmdArg MarkerColor(Color_t color)
RooCmdArg Slice(const RooArgSet &sliceSet)
const Bool_t kTRUE
Definition: Rtypes.h:91