Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf708_bphysics.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ##
5 ## Special p.d.f.'s: special decay pdf for B physics with mixing and/or CP violation
6 ##
7 ## \macro_code
8 ##
9 ## \date February 2018
10 ## \author Clemens Lange, Wouter Verkerke (C++ version)
11 
12 import ROOT
13 
14 
15 # B-decay with mixing
16 # -------------------------
17 
18 # Construct pdf
19 # -------------------------
20 
21 # Observable
22 dt = ROOT.RooRealVar("dt", "dt", -10, 10)
23 dt.setBins(40)
24 
25 # Parameters
26 dm = ROOT.RooRealVar("dm", "delta m(B0)", 0.472)
27 tau = ROOT.RooRealVar("tau", "tau (B0)", 1.547)
28 w = ROOT.RooRealVar("w", "flavour mistag rate", 0.1)
29 dw = ROOT.RooRealVar("dw", "delta mistag rate for B0/B0bar", 0.1)
30 
31 mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state")
32 mixState.defineType("mixed", -1)
33 mixState.defineType("unmixed", 1)
34 
35 tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0")
36 tagFlav.defineType("B0", 1)
37 tagFlav.defineType("B0bar", -1)
38 
39 # Use delta function resolution model
40 tm = ROOT.RooTruthModel("tm", "truth model", dt)
41 
42 # Construct Bdecay with mixing
43 bmix = ROOT.RooBMixDecay(
44  "bmix",
45  "decay",
46  dt,
47  mixState,
48  tagFlav,
49  tau,
50  dm,
51  w,
52  dw,
53  tm,
54  ROOT.RooBMixDecay.DoubleSided)
55 
56 # Plot pdf in various slices
57 # ---------------------------------------------------
58 
59 # Generate some data
60 data = bmix.generate(ROOT.RooArgSet(dt, mixState, tagFlav), 10000)
61 
62 # Plot B0 and B0bar tagged data separately
63 # For all plots below B0 and B0 tagged data will look somewhat differently
64 # if the flavor tagging mistag rate for B0 and B0 is different (i.e. dw!=0)
65 frame1 = dt.frame(ROOT.RooFit.Title(
66  "B decay distribution with mixing (B0/B0bar)"))
67 
68 data.plotOn(frame1, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
69 bmix.plotOn(frame1, ROOT.RooFit.Slice(tagFlav, "B0"))
70 
71 data.plotOn(frame1, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
72  ROOT.RooFit.MarkerColor(ROOT.kCyan))
73 bmix.plotOn(frame1, ROOT.RooFit.Slice(tagFlav, "B0bar"),
74  ROOT.RooFit.LineColor(ROOT.kCyan))
75 
76 # Plot mixed slice for B0 and B0bar tagged data separately
77 frame2 = dt.frame(ROOT.RooFit.Title(
78  "B decay distribution of mixed events (B0/B0bar)"))
79 
80 data.plotOn(frame2, ROOT.RooFit.Cut(
81  "mixState==mixState::mixed&&tagFlav==tagFlav::B0"))
82 bmix.plotOn(frame2, ROOT.RooFit.Slice(tagFlav, "B0"),
83  ROOT.RooFit.Slice(mixState, "mixed"))
84 
85 data.plotOn(
86  frame2,
87  ROOT.RooFit.Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),
88  ROOT.RooFit.MarkerColor(
89  ROOT.kCyan))
90 bmix.plotOn(frame2, ROOT.RooFit.Slice(tagFlav, "B0bar"), ROOT.RooFit.Slice(
91  mixState, "mixed"), ROOT.RooFit.LineColor(ROOT.kCyan))
92 
93 # Plot unmixed slice for B0 and B0bar tagged data separately
94 frame3 = dt.frame(ROOT.RooFit.Title(
95  "B decay distribution of unmixed events (B0/B0bar)"))
96 
97 data.plotOn(frame3, ROOT.RooFit.Cut(
98  "mixState==mixState::unmixed&&tagFlav==tagFlav::B0"))
99 bmix.plotOn(frame3, ROOT.RooFit.Slice(tagFlav, "B0"),
100  ROOT.RooFit.Slice(mixState, "unmixed"))
101 
102 data.plotOn(
103  frame3,
104  ROOT.RooFit.Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),
105  ROOT.RooFit.MarkerColor(
106  ROOT.kCyan))
107 bmix.plotOn(frame3, ROOT.RooFit.Slice(tagFlav, "B0bar"), ROOT.RooFit.Slice(
108  mixState, "unmixed"), ROOT.RooFit.LineColor(ROOT.kCyan))
109 
110 # B-decay with CP violation
111 # -------------------------
112 
113 # Construct pdf
114 # -------------------------
115 
116 # Additional parameters needed for B decay with CPV
117 CPeigen = ROOT.RooRealVar("CPeigen", "CP eigen value", -1)
118 absLambda = ROOT.RooRealVar("absLambda", "|lambda|", 1, 0, 2)
119 argLambda = ROOT.RooRealVar("absLambda", "|lambda|", 0.7, -1, 1)
120 effR = ROOT.RooRealVar("effR", "B0/B0bar reco efficiency ratio", 1)
121 
122 # Construct Bdecay with CP violation
123 bcp = ROOT.RooBCPEffDecay(
124  "bcp",
125  "bcp",
126  dt,
127  tagFlav,
128  tau,
129  dm,
130  w,
131  CPeigen,
132  absLambda,
133  argLambda,
134  effR,
135  dw,
136  tm,
137  ROOT.RooBCPEffDecay.DoubleSided)
138 
139 # Plot scenario 1 - sin(2b)=0.7, |l|=1
140 # ---------------------------------------------------------------------------
141 
142 # Generate some data
143 data2 = bcp.generate(ROOT.RooArgSet(dt, tagFlav), 10000)
144 
145 # Plot B0 and B0bar tagged data separately
146 frame4 = dt.frame(ROOT.RooFit.Title(
147  "B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)"))
148 
149 data2.plotOn(frame4, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
150 bcp.plotOn(frame4, ROOT.RooFit.Slice(tagFlav, "B0"))
151 
152 data2.plotOn(frame4, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
153  ROOT.RooFit.MarkerColor(ROOT.kCyan))
154 bcp.plotOn(frame4, ROOT.RooFit.Slice(tagFlav, "B0bar"),
155  ROOT.RooFit.LineColor(ROOT.kCyan))
156 
157 # # Plot scenario 2 - sin(2b)=0.7, |l|=0.7
158 # -------------------------------------------------------------------------------
159 
160 absLambda.setVal(0.7)
161 
162 # Generate some data
163 data3 = bcp.generate(ROOT.RooArgSet(dt, tagFlav), 10000)
164 
165 # Plot B0 and B0bar tagged data separately (sin2b = 0.7 plus direct CPV
166 # |l|=0.5)
167 frame5 = dt.frame(ROOT.RooFit.Title(
168  "B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)"))
169 
170 data3.plotOn(frame5, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
171 bcp.plotOn(frame5, ROOT.RooFit.Slice(tagFlav, "B0"))
172 
173 data3.plotOn(frame5, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
174  ROOT.RooFit.MarkerColor(ROOT.kCyan))
175 bcp.plotOn(frame5, ROOT.RooFit.Slice(tagFlav, "B0bar"),
176  ROOT.RooFit.LineColor(ROOT.kCyan))
177 
178 
179 # Generic B-decay with user coefficients
180 # -------------------------
181 
182 # Construct pdf
183 # -------------------------
184 
185 # Model parameters
186 DGbG = ROOT.RooRealVar("DGbG", "DGamma/GammaAvg", 0.5, -1, 1)
187 Adir = ROOT.RooRealVar("Adir", "-[1-abs(l)**2]/[1+abs(l)**2]", 0)
188 Amix = ROOT.RooRealVar("Amix", "2Im(l)/[1+abs(l)**2]", 0.7)
189 Adel = ROOT.RooRealVar("Adel", "2Re(l)/[1+abs(l)**2]", 0.7)
190 
191 # Derived input parameters for pdf
192 DG = ROOT.RooFormulaVar("DG", "Delta Gamma", "@1/@0",
193  ROOT.RooArgList(tau, DGbG))
194 
195 # Construct coefficient functions for sin,cos, modulations of decay
196 # distribution
197 fsin = ROOT.RooFormulaVar(
198  "fsin", "fsin", "@0*@1*(1-2*@2)", ROOT.RooArgList(Amix, tagFlav, w))
199 fcos = ROOT.RooFormulaVar(
200  "fcos", "fcos", "@0*@1*(1-2*@2)", ROOT.RooArgList(Adir, tagFlav, w))
201 fsinh = ROOT.RooFormulaVar("fsinh", "fsinh", "@0", ROOT.RooArgList(Adel))
202 
203 # Construct generic B decay pdf using above user coefficients
204 bcpg = ROOT.RooBDecay("bcpg", "bcpg", dt, tau, DG, ROOT.RooFit.RooConst(
205  1), fsinh, fcos, fsin, dm, tm, ROOT.RooBDecay.DoubleSided)
206 
207 # Plot - Im(l)=0.7, e(l)=0.7 |l|=1, G/G=0.5
208 # -------------------------------------------------------------------------------------
209 
210 # Generate some data
211 data4 = bcpg.generate(ROOT.RooArgSet(dt, tagFlav), 10000)
212 
213 # Plot B0 and B0bar tagged data separately
214 frame6 = dt.frame(ROOT.RooFit.Title(
215  "B decay distribution with CPV(Im(l)=0.7,Re(l)=0.7,|l|=1,dG/G=0.5) (B0/B0bar)"))
216 
217 data4.plotOn(frame6, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
218 bcpg.plotOn(frame6, ROOT.RooFit.Slice(tagFlav, "B0"))
219 
220 data4.plotOn(frame6, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
221  ROOT.RooFit.MarkerColor(ROOT.kCyan))
222 bcpg.plotOn(frame6, ROOT.RooFit.Slice(tagFlav, "B0bar"),
223  ROOT.RooFit.LineColor(ROOT.kCyan))
224 
225 c = ROOT.TCanvas("rf708_bphysics", "rf708_bphysics", 1200, 800)
226 c.Divide(3, 2)
227 c.cd(1)
228 ROOT.gPad.SetLeftMargin(0.15)
229 frame1.GetYaxis().SetTitleOffset(1.6)
230 frame1.Draw()
231 c.cd(2)
232 ROOT.gPad.SetLeftMargin(0.15)
233 frame2.GetYaxis().SetTitleOffset(1.6)
234 frame2.Draw()
235 c.cd(3)
236 ROOT.gPad.SetLeftMargin(0.15)
237 frame3.GetYaxis().SetTitleOffset(1.6)
238 frame3.Draw()
239 c.cd(4)
240 ROOT.gPad.SetLeftMargin(0.15)
241 frame4.GetYaxis().SetTitleOffset(1.6)
242 frame4.Draw()
243 c.cd(5)
244 ROOT.gPad.SetLeftMargin(0.15)
245 frame5.GetYaxis().SetTitleOffset(1.6)
246 frame5.Draw()
247 c.cd(6)
248 ROOT.gPad.SetLeftMargin(0.15)
249 frame6.GetYaxis().SetTitleOffset(1.6)
250 frame6.Draw()
251 
252 c.SaveAs("rf708_bphysics.png")