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