Logo ROOT  
Reference Guide
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## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14
15# B-decay with mixing
16# -------------------------
17
18# Construct pdf
19# -------------------------
20
21# Observable
22dt = ROOT.RooRealVar("dt", "dt", -10, 10)
23dt.setBins(40)
24
25# Parameters
26dm = ROOT.RooRealVar("dm", "delta m(B0)", 0.472)
27tau = ROOT.RooRealVar("tau", "tau (B0)", 1.547)
28w = ROOT.RooRealVar("w", "flavour mistag rate", 0.1)
29dw = ROOT.RooRealVar("dw", "delta mistag rate for B0/B0bar", 0.1)
30
31mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state")
32mixState.defineType("mixed", -1)
33mixState.defineType("unmixed", 1)
34
35tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0")
36tagFlav.defineType("B0", 1)
37tagFlav.defineType("B0bar", -1)
38
39# Use delta function resolution model
40tm = ROOT.RooTruthModel("tm", "truth model", dt)
41
42# Construct Bdecay with mixing
43bmix = 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
60data = 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)
65frame1 = dt.frame(ROOT.RooFit.Title(
66 "B decay distribution with mixing (B0/B0bar)"))
67
68data.plotOn(frame1, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
69bmix.plotOn(frame1, ROOT.RooFit.Slice(tagFlav, "B0"))
70
71data.plotOn(frame1, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
72 ROOT.RooFit.MarkerColor(ROOT.kCyan))
73bmix.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
77frame2 = dt.frame(ROOT.RooFit.Title(
78 "B decay distribution of mixed events (B0/B0bar)"))
79
80data.plotOn(frame2, ROOT.RooFit.Cut(
81 "mixState==mixState::mixed&&tagFlav==tagFlav::B0"))
82bmix.plotOn(frame2, ROOT.RooFit.Slice(tagFlav, "B0"),
83 ROOT.RooFit.Slice(mixState, "mixed"))
84
85data.plotOn(
86 frame2,
87 ROOT.RooFit.Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),
88 ROOT.RooFit.MarkerColor(
89 ROOT.kCyan))
90bmix.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
94frame3 = dt.frame(ROOT.RooFit.Title(
95 "B decay distribution of unmixed events (B0/B0bar)"))
96
97data.plotOn(frame3, ROOT.RooFit.Cut(
98 "mixState==mixState::unmixed&&tagFlav==tagFlav::B0"))
99bmix.plotOn(frame3, ROOT.RooFit.Slice(tagFlav, "B0"),
100 ROOT.RooFit.Slice(mixState, "unmixed"))
101
102data.plotOn(
103 frame3,
104 ROOT.RooFit.Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),
105 ROOT.RooFit.MarkerColor(
106 ROOT.kCyan))
107bmix.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
117CPeigen = ROOT.RooRealVar("CPeigen", "CP eigen value", -1)
118absLambda = ROOT.RooRealVar("absLambda", "|lambda|", 1, 0, 2)
119argLambda = ROOT.RooRealVar("absLambda", "|lambda|", 0.7, -1, 1)
120effR = ROOT.RooRealVar("effR", "B0/B0bar reco efficiency ratio", 1)
121
122# Construct Bdecay with CP violation
123bcp = 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
143data2 = bcp.generate(ROOT.RooArgSet(dt, tagFlav), 10000)
144
145# Plot B0 and B0bar tagged data separately
146frame4 = dt.frame(ROOT.RooFit.Title(
147 "B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)"))
148
149data2.plotOn(frame4, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
150bcp.plotOn(frame4, ROOT.RooFit.Slice(tagFlav, "B0"))
151
152data2.plotOn(frame4, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
153 ROOT.RooFit.MarkerColor(ROOT.kCyan))
154bcp.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
160absLambda.setVal(0.7)
161
162# Generate some data
163data3 = 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)
167frame5 = dt.frame(ROOT.RooFit.Title(
168 "B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)"))
169
170data3.plotOn(frame5, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
171bcp.plotOn(frame5, ROOT.RooFit.Slice(tagFlav, "B0"))
172
173data3.plotOn(frame5, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
174 ROOT.RooFit.MarkerColor(ROOT.kCyan))
175bcp.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
186DGbG = ROOT.RooRealVar("DGbG", "DGamma/GammaAvg", 0.5, -1, 1)
187Adir = ROOT.RooRealVar("Adir", "-[1-abs(l)**2]/[1+abs(l)**2]", 0)
188Amix = ROOT.RooRealVar("Amix", "2Im(l)/[1+abs(l)**2]", 0.7)
189Adel = ROOT.RooRealVar("Adel", "2Re(l)/[1+abs(l)**2]", 0.7)
190
191# Derived input parameters for pdf
192DG = 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
197fsin = ROOT.RooFormulaVar(
198 "fsin", "fsin", "@0*@1*(1-2*@2)", ROOT.RooArgList(Amix, tagFlav, w))
199fcos = ROOT.RooFormulaVar(
200 "fcos", "fcos", "@0*@1*(1-2*@2)", ROOT.RooArgList(Adir, tagFlav, w))
201fsinh = ROOT.RooFormulaVar("fsinh", "fsinh", "@0", ROOT.RooArgList(Adel))
202
203# Construct generic B decay pdf using above user coefficients
204bcpg = 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
211data4 = bcpg.generate(ROOT.RooArgSet(dt, tagFlav), 10000)
212
213# Plot B0 and B0bar tagged data separately
214frame6 = 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
217data4.plotOn(frame6, ROOT.RooFit.Cut("tagFlav==tagFlav::B0"))
218bcpg.plotOn(frame6, ROOT.RooFit.Slice(tagFlav, "B0"))
219
220data4.plotOn(frame6, ROOT.RooFit.Cut("tagFlav==tagFlav::B0bar"),
221 ROOT.RooFit.MarkerColor(ROOT.kCyan))
222bcpg.plotOn(frame6, ROOT.RooFit.Slice(tagFlav, "B0bar"),
223 ROOT.RooFit.LineColor(ROOT.kCyan))
224
225c = ROOT.TCanvas("rf708_bphysics", "rf708_bphysics", 1200, 800)
226c.Divide(3, 2)
227c.cd(1)
228ROOT.gPad.SetLeftMargin(0.15)
229frame1.GetYaxis().SetTitleOffset(1.6)
230frame1.Draw()
231c.cd(2)
232ROOT.gPad.SetLeftMargin(0.15)
233frame2.GetYaxis().SetTitleOffset(1.6)
234frame2.Draw()
235c.cd(3)
236ROOT.gPad.SetLeftMargin(0.15)
237frame3.GetYaxis().SetTitleOffset(1.6)
238frame3.Draw()
239c.cd(4)
240ROOT.gPad.SetLeftMargin(0.15)
241frame4.GetYaxis().SetTitleOffset(1.6)
242frame4.Draw()
243c.cd(5)
244ROOT.gPad.SetLeftMargin(0.15)
245frame5.GetYaxis().SetTitleOffset(1.6)
246frame5.Draw()
247c.cd(6)
248ROOT.gPad.SetLeftMargin(0.15)
249frame6.GetYaxis().SetTitleOffset(1.6)
250frame6.Draw()
251
252c.SaveAs("rf708_bphysics.png")