Belle II Software  release-05-01-25
B0_GenDeltaTFit.py
1 #!/usr/bin/env python
2 # -*- coding: utf-8 -*-
3 
4 
21 
22 import ROOT
23 import basf2 as b2
24 from ROOT import Belle2
25 import modularAnalysis as ma
26 import numpy as np
27 import matplotlib.pyplot as plt
28 import pylab
29 import sys
30 import math
31 import random
32 from variables import variables as mc_variables
33 
34 inputFiles = []
35 Xlim1 = 0
36 Xlim2 = 999
37 
38 validation_path = b2.Path()
39 
40 ma.inputMdstList("default", inputFiles, path=validation_path)
41 
42 B0s = ('B0:all', '')
43 ma.fillParticleListsFromMC([B0s], True, path=validation_path)
44 
45 ma.matchMCTruth('B0:all', path=validation_path)
46 
47 Upsilon4S = ('Upsilon(4S):all', '')
48 ma.fillParticleListsFromMC([Upsilon4S], True, path=validation_path)
49 
50 ma.matchMCTruth('Upsilon(4S):all', path=validation_path)
51 
52 labelAsymmetry = "#frac{N^{#it{B}^{0}} - N^{#bar{#it{B}}^{0}}}{N^{#it{B}^{0}} + N^{#bar{#it{B}}^{0}}}"
53 
54 
55 def plotWithAsymmetry(
56  rooFitFrame,
57  rooRealVar,
58  dotsAll,
59  dotsA,
60  dotsB,
61  curveAll,
62  curveA,
63  curveB,
64  units,
65  nameOfPlot,
66  legend,
67  textFitRes,
68  asymmVarName,
69  removeArtifacts):
70  """
71  Produces a plot with an asymmetry (y_A -y_B)/(y_A + y_B) below, where y_A and y_B are the y values of curves A and B.
72  """
73 
74  # plot for Thesis
75  rooFitFrame.Print()
76  rooFitFrame.GetXaxis().SetTitle("")
77  rooFitFrame.GetXaxis().SetLabelSize(0)
78 
79  rooFitFrame.GetYaxis().SetTitleSize(0.072)
80  rooFitFrame.GetYaxis().SetTitleOffset(0.98)
81  rooFitFrame.GetYaxis().SetLabelSize(0.055)
82 
83  # yLabelBin = 0
84  # yLabelBin = '{:0.0f}'.format(float(dots.GetErrorXlow(1) + dots.GetErrorXhigh(1)))
85 
86  pointsHist = ROOT.RooHist()
87 
88  iBin = 0
89 
90  xValModel = ROOT.Double(-1.E30)
91  yValModel = ROOT.Double(-1.E30)
92  xValDot = ROOT.Double(-1.E30)
93  yValDotA = ROOT.Double(-1.E30)
94  yValDotB = ROOT.Double(-1.E30)
95 
96  iDotsA = ROOT.RooRealVar("iDotsA", "", 0.)
97  idotsB = ROOT.RooRealVar("idotsB", "", 0.)
98  iDotError = ROOT.RooRealVar("iDotErrorA", "", 0.)
99  # yComb = ROOT.RooFormulaVar("yComb", "", "-1*@0", ROOT.RooArgList(y1RV))
100  iAsymmetry = ROOT.RooFormulaVar("yComb", "", "(@0 - @1)/(@0 + @1)", ROOT.RooArgList(iDotsA, idotsB))
101 
102  while iBin < dotsA.GetN() - 1:
103  dotsA.GetPoint(iBin, xValDot, yValDotA)
104  dotsB.GetPoint(iBin, xValDot, yValDotB)
105 
106  iDotsA.setVal(yValDotA)
107  idotsB.setVal(yValDotB)
108 
109  errA = float(dotsA.GetErrorYlow(iBin) + dotsA.GetErrorYhigh(iBin))
110  errB = float(dotsB.GetErrorYlow(iBin) + dotsB.GetErrorYhigh(iBin))
111 
112  if np.isnan(iAsymmetry.getVal()) is not True:
113 
114  asymValue = iAsymmetry.getVal()
115 
116  if (yValDotA + yValDotB) > 0:
117  asymError = float(2 * math.sqrt(yValDotA**2 * errA**2 + yValDotB**2 * errB**2) / (yValDotA + yValDotB)**2)
118  else:
119  asymError = 0
120 
121  if removeArtifacts:
122  if abs(asymValue) > 3:
123  asymValue = 0
124 
125  # print("xval = ", xValDot)
126  # print(iBin, " = ", iAsymmetry.getVal(), ", @0 = ", yValDot, ", @1 = ", curveA.interpolate(xValDot), ", @2 = ",
127  # float(dots.GetErrorYlow(iBin) + dots.GetErrorYhigh(iBin)))
128  pointsHist.addBinWithXYError(xValDot, asymValue,
129  dotsA.GetErrorXlow(iBin),
130  dotsA.GetErrorXhigh(iBin),
131  asymError, asymError,)
132 
133  iBin += 1
134 
135  pointsHist.SetMarkerStyle(dotsA.GetMarkerStyle())
136  # pointsHist.SetMarkerSize(dotsA.GetMarkerSize())
137 
138  rooFitFramAsym = rooRealVar.frame()
139  # rooFitFramAsym.addObject(mbcFrame.pullHist())
140  # rooFitFramAsym.addPlotable(pointsHist, "P")
141 
142  rooFitFramAsym.GetXaxis().SetTickSize(0.07)
143  rooFitFramAsym.GetYaxis().SetTickSize(0.024)
144 
145  rooFitFramAsym.SetTitle("")
146  rooFitFramAsym.GetXaxis().SetTitle(rooRealVar.GetTitle() + " " + units)
147  rooFitFramAsym.GetYaxis().SetTitle(asymmVarName + " ")
148 
149  if asymmVarName == "#it{a}(#Delta#it{t})":
150  rooFitFramAsym.GetXaxis().SetTitleSize(0.2)
151  rooFitFramAsym.GetXaxis().SetTitleOffset(0.9)
152  rooFitFramAsym.GetYaxis().SetTitleSize(0.19)
153  rooFitFramAsym.GetYaxis().SetTitleOffset(0.32)
154  else:
155  rooFitFramAsym.GetXaxis().SetTitleSize(0.18)
156  rooFitFramAsym.GetXaxis().SetTitleOffset(0.95)
157  rooFitFramAsym.GetYaxis().SetTitleSize(0.13)
158  rooFitFramAsym.GetYaxis().SetTitleOffset(0.52)
159 
160  rooFitFramAsym.GetXaxis().SetLabelSize(0.155)
161  rooFitFramAsym.GetYaxis().SetLabelSize(0.120)
162 
163  rooFitFramAsym.SetAxisRange(-1.0, 1.0, "Y")
164  rooFitFramAsym.GetYaxis().SetNdivisions(10)
165  rooFitFramAsym.GetYaxis().ChangeLabel(1, -1, 0.)
166  rooFitFramAsym.GetYaxis().ChangeLabel(3, -1, 0.)
167  rooFitFramAsym.GetYaxis().ChangeLabel(5, -1, 0.)
168  rooFitFramAsym.GetYaxis().ChangeLabel(6, -1, -1, -1, -1, -1, "0.0")
169  rooFitFramAsym.GetYaxis().ChangeLabel(7, -1, 0.)
170  rooFitFramAsym.GetYaxis().ChangeLabel(9, -1, 0.)
171  rooFitFramAsym.GetYaxis().ChangeLabel(11, -1, 0.)
172 
173  # print("Number of Y axis bins = ", rooFitFramAsym.GetYaxis().GetNbins() )
174 
175  gLine0 = ROOT.TLine(rooFitFramAsym.GetXaxis().GetXmin(), 0,
176  rooFitFramAsym.GetXaxis().GetXmax(), 0)
177  gLine0.SetLineWidth(1)
178  gLine0.SetLineColor(13)
179  gLine0.SetLineStyle(3)
180 
181  gLine1 = ROOT.TLine(rooFitFramAsym.GetXaxis().GetXmin(), 0.4,
182  rooFitFramAsym.GetXaxis().GetXmax(), 0.4)
183  gLine1.SetLineWidth(1)
184  gLine1.SetLineColor(13)
185  gLine1.SetLineStyle(3)
186 
187  gLine2 = ROOT.TLine(rooFitFramAsym.GetXaxis().GetXmin(), -0.4,
188  rooFitFramAsym.GetXaxis().GetXmax(), -0.4)
189  gLine2.SetLineWidth(1)
190  gLine2.SetLineColor(13)
191  gLine2.SetLineStyle(3)
192 
193  gLine3 = ROOT.TLine(rooFitFramAsym.GetXaxis().GetXmin(), 0.8,
194  rooFitFramAsym.GetXaxis().GetXmax(), 0.8)
195  gLine3.SetLineWidth(1)
196  gLine3.SetLineColor(13)
197  gLine3.SetLineStyle(3)
198 
199  gLine4 = ROOT.TLine(rooFitFramAsym.GetXaxis().GetXmin(), -0.8,
200  rooFitFramAsym.GetXaxis().GetXmax(), -0.8)
201  gLine4.SetLineWidth(1)
202  gLine4.SetLineColor(13)
203  gLine4.SetLineStyle(3)
204 
205  c1 = ROOT.TCanvas("c1", "c1", 700, 640)
206  c1.SetBottomMargin(0)
207  c1.cd()
208 
209  Pad1 = ROOT.TPad("p1", "p1", 0, 0.277, 1, 1, 0)
210  Pad2 = ROOT.TPad("p2", "p2", 0, 0, 1, 0.276, 0)
211  Pad1.Draw()
212  Pad2.Draw()
213 
214  Pad1.SetLeftMargin(0.15)
215  Pad1.SetBottomMargin(0.02)
216  Pad1.SetTopMargin(0.06)
217 
218  Pad2.SetLeftMargin(0.15)
219  Pad2.SetBottomMargin(0.4)
220  Pad2.SetTopMargin(0.01)
221 
222  Pad2.cd()
223  rooFitFramAsym.Draw()
224  gLine0.Draw("SAME")
225  gLine1.Draw("SAME")
226  gLine2.Draw("SAME")
227  gLine3.Draw("SAME")
228  gLine4.Draw("SAME")
229  pointsHist.Draw("P SAME")
230  # Pad2.Update()
231 
232  Pad1.cd()
233 
234  rooFitFrame.Draw()
235 
236  # if textFitRes != "":
237  #
238  # lFitRes = ROOT.TLegend(0.238, 0.95, 0.7, 0.98)
239  # lFitRes.AddEntry("", textFitRes , " ")
240  # lFitRes.SetBorderSize(0)
241  # lFitRes.SetTextSize(0.056)
242  # lFitRes.Draw()
243 
244  if legend != "":
245  legend.Draw()
246  # Pad1.Update()
247  nameOfPlot = "WithAsymm" + nameOfPlot + ".pdf"
248  c1.SaveAs(nameOfPlot)
249  c1.Clear()
250 
251 
252 class fitDeltaT(b2.Module):
253  """
254  This class makes plots to validate the generated CP asymmetry. It saves the decay time information of the simulation
255  for both B mesons by running as a basf2 module. It fills the information in RooDataSets and at the end it fits the
256  known quantum mechanical pdfs to the delta T distribution and the individual decay time distributions of
257  the signal and the tag-side B mesons to check if the S and the A parameters in the simulation are correct.
258  """
259 
323 
324  def __init__(self):
325  """
326  Here the module initializes by declaring the variables corresponding to the interesting physical parameters and
327  decay times. The RooDataSets where the information will be saved are also declared.
328  """
329  super(fitDeltaT, self).__init__()
330  self.nentries = 0
331 
332  self.DT = ROOT.RooRealVar("DT", "#Delta#it{t}", 0., -11., 11., "ps")
333  self.Tsig = ROOT.RooRealVar("Tsig", "#it{t}_{sig}^{gen}", 0., 0., 11., "ps")
334  self.Ttag = ROOT.RooRealVar("Ttag", "#it{t}_{tag}^{gen}", 0., 0., 11., "ps")
335  self.DT0 = ROOT.RooRealVar("DT0", "DT0", 0., -7., 7.)
336  self.Tau = ROOT.RooRealVar("Tau", "Tau", 1.525, 0., 4.)
337  self.A = ROOT.RooRealVar("A", "A", 0, -1, 1)
338  self.S = ROOT.RooRealVar("S", "S", 0, -1, 1)
339  self.DM = ROOT.RooRealVar("DM", "DM", 0.507)
340  self.Norm = ROOT.RooRealVar("Norm", "Norm", 2, 1, 8)
341  self.BetaGamma = ROOT.RooRealVar("BetaGamma", "BetaGamma", 0.5, 0.1, 1)
342 
343  self.q = ROOT.RooCategory("q", "q")
344  self.q.defineType("1")
345  self.q.defineType("-1")
346 
347  self.qsig = ROOT.RooCategory("qsig", "qsig")
348  self.qsig.defineType("1")
349  self.qsig.defineType("-1")
350 
351  self.fitData = ROOT.RooDataSet("fitData", "fitData", ROOT.RooArgSet(self.DT, self.q))
352  self.fitDataTag = ROOT.RooDataSet("fitDataTag", "fitDataTag", ROOT.RooArgSet(self.Ttag, self.q))
353  self.fitDataTagDeltaTPos = ROOT.RooDataSet("fitDataTagDeltaTPos", "fitDataTagDeltaTPos", ROOT.RooArgSet(self.Ttag, self.q))
354  self.fitDataTagDeltaTNeg = ROOT.RooDataSet("fitDataTagDeltaTNeg", "fitDataTagDeltaTNeg", ROOT.RooArgSet(self.Ttag, self.q))
355  self.fitDataSig = ROOT.RooDataSet("fitDataSig", "fitDataSig", ROOT.RooArgSet(self.Tsig, self.q, self.qsig))
356  self.fitDataSigDeltaTPos = ROOT.RooDataSet("fitDataSigDeltaTPos", "fitDataSigDeltaTPos", ROOT.RooArgSet(self.Tsig, self.q))
357  self.fitDataSigDeltaTNeg = ROOT.RooDataSet("fitDataSigDeltaTNeg", "fitDataSigDeltaTNeg", ROOT.RooArgSet(self.Tsig, self.q))
358 
359  self.B0sInTagSide = 0
361 
364  self.fractionB0barInSignalSide = 0.7057091319609885
365  self.fractionB0InSignalSide = 0.2942908680390115
366 
367  # Study od Delta z
368 
369  self.DeltaZsigUpsilon = ROOT.RooRealVar("DeltaZsigUpsilon", "(#it{z}_{sig}^{dec} - #it{z}_{sig}^{prod})^{gen}", 0., 0., 0.8)
370  self.DeltaZtagUpsilon = ROOT.RooRealVar(
371  "DeltaZtagUpsilon",
372  "(#it{z}_{tag}^{dec} - #it{z}_{tag}^{prod})^{gen}",
373  0.,
374  0.,
375  0.8,
376  "mm")
377  self.fitDataDeltaZsigUpsilon = ROOT.RooDataSet(
378  "fitDataDeltaZsigUpsilon",
379  "fitDataDeltaZsigUpsilon",
380  ROOT.RooArgSet(
381  self.DeltaZsigUpsilon,
382  self.q))
383  self.fitDataDeltaZtagUpsilon = ROOT.RooDataSet(
384  "fitDataDeltaZtagUpsilon",
385  "fitDataDeltaZtagUpsilon",
386  ROOT.RooArgSet(
387  self.DeltaZtagUpsilon,
388  self.q))
389 
390  # Only for plotting
391 
392  self.TsigPosNeg = ROOT.RooRealVar("Tsig", "#it{t}_{sig}^{gen}", 0., 0., 7., "ps")
393  self.TtagPosNeg = ROOT.RooRealVar("Ttag", "#it{t}_{tag}^{gen}", 0., 0., 7., "ps")
394 
395  def event(self):
396  """
397  Here the information is collected event by event.
398  """
399 
400  plist = Belle2.PyStoreObj('B0:all')
401  plistUpsilon = Belle2.PyStoreObj('Upsilon(4S):all')
402 
403  tagPDG = 0
404  sigPDG = 0
405  tSig = 0
406  tTag = 0
407 
408  zSig = 0
409  zTag = 0
410  zUpsilon = 0
411 
412  DeltaT = 0
413  B0Flags = []
414 
415  for index in range(0, plistUpsilon.obj().getListSize()):
416  Upsilon4Sparticle = plistUpsilon.obj().getParticle(index)
417 
418  if Upsilon4Sparticle.getPDGCode() == 300553:
419  zUpsilon = mc_variables.evaluate("mcZ", Upsilon4Sparticle)
420 
421  for index in range(0, plist.obj().getListSize()):
422  B0particle = plist.obj().getParticle(index) # Pointer to the particle with highest prob
423  # MCB0particle = B0particle.getMCParticle()
424 
425  isB0sig = False
426 
427  if B0particle.getNDaughters() == 2:
428  B0daughters = []
429  for B0daughter in B0particle.getDaughters():
430  B0daughters.append(B0daughter)
431  if (abs(B0daughters[0].getPDGCode()) == 443 and abs(B0daughters[1].getPDGCode()) == 310) or \
432  (abs(B0daughters[1].getPDGCode()) == 443 and abs(B0daughters[0].getPDGCode()) == 310) or \
433  (abs(B0daughters[1].getPDGCode()) == 16 and abs(B0daughters[0].getPDGCode()) == 16):
434  isB0sig = True
435  # print("hi")
436 
437  if isB0sig:
438  tSig = mc_variables.evaluate("mcDecayTime", B0particle)
439  zSig = mc_variables.evaluate("mcZ", B0particle)
440  sigPDG = B0particle.getPDGCode()
441  else:
442  tTag = mc_variables.evaluate("mcDecayTime", B0particle)
443  zTag = mc_variables.evaluate("mcZ", B0particle)
444  tagPDG = B0particle.getPDGCode()
445 
446  # print(tSig, " ", tTag)
447  B0Flags.append(isB0sig)
448 
449  tagProb = random.uniform(0, 1)
450  signalProb = random.uniform(0, 1)
451 
452  rejectEventFlag = False
453 
454  DeltaT = 1000 * (tSig - tTag)
455  # DeltaZ = (zSig - zTag)
456  # DeltaT = 100*DeltaZ /(0.284*3)
457  if len(B0Flags) == 2:
458  if B0Flags[0] != B0Flags[1] and rejectEventFlag is not True:
459  self.nentries += 1
460  # DeltaT = 1000*(tSig - tTag)
461  # print(DeltaT)
462  self.DT.setVal(DeltaT)
463  self.Tsig.setVal(float(1000 * tSig))
464  self.Ttag.setVal(float(1000 * tTag))
465 
466  if tagPDG > 0:
467  self.q.setLabel("1")
468  self.B0sInTagSide += 1
469  elif tagPDG < 0:
470  self.q.setLabel("-1")
471  self.B0barsInTagSide += 1
472  if sigPDG > 0:
473  self.qsig.setLabel("1")
474  self.B0sInSignalSide += 1
475  elif sigPDG < 0:
476  self.qsig.setLabel("-1")
477  self.B0barsInSignalSide += 1
478  self.fitData.add(ROOT.RooArgSet(self.DT, self.q))
479  self.fitDataSig.add(ROOT.RooArgSet(self.Tsig, self.q, self.qsig))
480  self.fitDataTag.add(ROOT.RooArgSet(self.Ttag, self.q))
481 
482  if DeltaT >= 0:
483  self.fitDataTagDeltaTPos.add(ROOT.RooArgSet(self.Ttag, self.q))
484  self.fitDataSigDeltaTPos.add(ROOT.RooArgSet(self.Tsig, self.q))
485 
486  if DeltaT <= 0:
487  self.fitDataTagDeltaTNeg.add(ROOT.RooArgSet(self.Ttag, self.q))
488  self.fitDataSigDeltaTNeg.add(ROOT.RooArgSet(self.Tsig, self.q))
489 
490  self.DeltaZtagUpsilon.setVal(10 * float(zTag - zUpsilon))
491  self.DeltaZsigUpsilon.setVal(10 * float(zSig - zUpsilon))
492  self.fitDataDeltaZtagUpsilon.add(ROOT.RooArgSet(self.DeltaZtagUpsilon, self.q))
493  self.fitDataDeltaZsigUpsilon.add(ROOT.RooArgSet(self.DeltaZsigUpsilon, self.q))
494 
495  def terminate(self):
496  """
497  Here the known quantum mechanical pdfs are defined and fitted to the generated MC distributions.
498  Afterwards, the plots are saved.
499  """
500  # S
501 
502  # Delta t Fit
503 
504  self.fitData.Print()
505 
506  self.DM.setConstant(ROOT.kTRUE)
507  self.Tau.setConstant(ROOT.kTRUE)
508  self.DT0.setConstant(ROOT.kTRUE)
509  self.Norm.setConstant(ROOT.kTRUE)
510  # self.A.setConstant(ROOT.kTRUE)
511  # self.S.setConstant(ROOT.kTRUE)
512 
513  # self.S.setVal(6.9184e-01)
514  # self.S.setConstant(ROOT.kTRUE)
515  # self.A.setVal(4.4702e-03)
516  # self.A.setConstant(ROOT.kTRUE)
517 
518  DeltaTp = ROOT.RooGenericPdf(
519  "DeltaTp",
520  "DeltaTp",
521  "(exp(-1*abs(DT-DT0)/Tau)/(2*Norm*Tau))*(1+(A*cos(DM*DT)+S*sin(DM*DT)))",
522  ROOT.RooArgList(
523  self.Norm,
524  self.DT,
525  self.DT0,
526  self.Tau,
527  self.A,
528  self.S,
529  self.DM))
530  DeltaTm = ROOT.RooGenericPdf(
531  "DeltaTm",
532  "DeltaTm",
533  "(exp(-1*abs(DT-DT0)/Tau)/(2*Norm*Tau))*(1-(A*cos(DM*DT)+S*sin(DM*DT)))",
534  ROOT.RooArgList(
535  self.Norm,
536  self.DT,
537  self.DT0,
538  self.Tau,
539  self.A,
540  self.S,
541  self.DM))
542  DeltaT = ROOT.RooSimultaneous("DeltaT", "DeltaT", self.q)
543  DeltaT.addPdf(DeltaTp, "1")
544  DeltaT.addPdf(DeltaTm, "-1")
545  DeltaT.Print()
546 
547  fitRes = DeltaT.fitTo(
548  self.fitData, ROOT.RooFit.Minos(
549  ROOT.kFALSE), ROOT.RooFit.Extended(
550  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
551  fitRes.Print("v")
552 
553  DTFrame = self.DT.frame()
554  print("here", str(DTFrame.GetXaxis().GetXmax()))
555  self.fitData.plotOn(DTFrame, ROOT.RooFit.Cut("abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())))
556  self.fitData.plotOn(DTFrame,
557  ROOT.RooFit.Cut("q==q::1" + " && abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())),
558  ROOT.RooFit.MarkerColor(ROOT.kRed),
559  ROOT.RooFit.Name("data_histB0"))
560  self.fitData.plotOn(DTFrame, ROOT.RooFit.Cut("q==1" + " && abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())),
561  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1), ROOT.RooFit.Name("data_histAB0"))
562 
563  # DeltaT.plotOn(DTFrame, ROOT.RooFit.Slice(self.q,"1"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q),self.fitData),
564  # ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor( ROOT.kRed),ROOT.RooFit.Components("B02JpsiKsp"));
565  # DeltaT.plotOn(DTFrame, ROOT.RooFit.Slice(self.q,"-1"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q),self.fitData),
566  # ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor( ROOT.kBlue+1), ROOT.RooFit.Components("B02JpsiKsm"));
567  DeltaT.plotOn(
568  DTFrame,
569  ROOT.RooFit.ProjWData(
570  ROOT.RooArgSet(
571  self.q),
572  self.fitData),
573  ROOT.RooFit.LineWidth(3),
574  ROOT.RooFit.LineColor(
575  ROOT.kBlack),
576  ROOT.RooFit.LineStyle(3))
577  DeltaT.plotOn(
578  DTFrame,
579  ROOT.RooFit.ProjWData(
580  ROOT.RooArgSet(
581  self.q),
582  self.fitData),
583  ROOT.RooFit.LineWidth(3),
584  ROOT.RooFit.LineColor(
585  ROOT.kRed),
586  ROOT.RooFit.Components("DeltaTp"))
587  DeltaT.plotOn(DTFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q), self.fitData), ROOT.RooFit.LineWidth(
588  3), ROOT.RooFit.LineColor(ROOT.kBlue + 1), ROOT.RooFit.LineStyle(7), ROOT.RooFit.Components("DeltaTm"))
589 
590  DTFrame.SetTitle("")
591  DTFrame.GetXaxis().SetTitle("#Delta#it{t}^{gen} [ps]")
592  DTFrame.GetXaxis().SetTitleSize(0.05)
593  DTFrame.GetXaxis().SetLabelSize(0.045)
594  DTFrame.GetYaxis().SetTitleSize(0.05)
595  DTFrame.GetYaxis().SetTitleOffset(1.5)
596  DTFrame.GetYaxis().SetLabelSize(0.045)
597 
598  c1 = ROOT.TCanvas("c1", "c1", 1400, 1100)
599  c1.cd()
600 
601  self.DT.Print()
602 
603  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
604  Pad.SetLeftMargin(0.15)
605  Pad.SetBottomMargin(0.15)
606  Pad.Draw()
607  Pad.cd()
608  DTFrame.Draw()
609  DTFrame.Print()
610 
611  curveMC = DTFrame.getCurve("DeltaT_Norm[DT]")
612  curveMCB0 = DTFrame.getCurve("DeltaT_Norm[DT]_Comp[DeltaTp]")
613  curveMCAB0 = DTFrame.getCurve("DeltaT_Norm[DT]_Comp[DeltaTm]")
614  curveMC.SetFillColor(ROOT.kWhite)
615  curveMCB0.SetFillColor(ROOT.kWhite)
616  curveMCAB0.SetFillColor(ROOT.kWhite)
617 
618  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
619  lMC.AddEntry(curveMC, "Both")
620  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
621  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
622 
623  lMC.SetTextSize(0.065)
624  lMC.Draw()
625 
626  textFitRes = "#it{A}_{#it{CP}} = " + \
627  '{: 1.2f}'.format(self.A.getValV() if abs(self.A.getValV()) > 0.005 else abs(self.A.getValV())) + \
628  ", " + "#it{S}_{#it{CP}} = " + \
629  '{: 1.2f}'.format(self.S.getValV() if abs(self.S.getValV()) > 0.005 else abs(self.S.getValV()))
630 
631  lFitRes = ROOT.TLegend(0.215, 0.905, 0.7, 0.98)
632  lFitRes.AddEntry("", textFitRes, " ")
633  lFitRes.SetBorderSize(0)
634  lFitRes.SetTextSize(0.054)
635  lFitRes.Draw()
636 
637  # lFitRes = ROOT.TLegend(0.18, 0.68, 0.40, 0.89)
638  # lFitRes.AddEntry("", "#it{A}_{#it{CP}} = " + '{: 1.2f}'.format(self.A.getValV()), " ")
639  # lFitRes.AddEntry("", "#it{S}_{#it{CP}} = " + '{: 1.2f}'.format(self.S.getValV()), " " )
640  # lFitRes.SetBorderSize(0)
641  # lFitRes.SetTextSize(0.054)
642  # lFitRes.Draw()
643 
644  textValA = ROOT.TText(
645  0.1,
646  1,
647  "#it{A}_{#it{CP}} = " +
648  '{:1.2f}'.format(
649  self.A.getValV()) +
650  " #pm " +
651  '{:1.2f}'.format(
652  self.A.getError()))
653  textValS = ROOT.TText(
654  0.1,
655  0.5,
656  "#it{A}_{#it{CP}} = " +
657  '{:1.2f}'.format(
658  self.S.getValV()) +
659  " #pm " +
660  '{:1.2f}'.format(
661  self.S.getError()))
662 
663  fitRes.Print("v")
664  c1.SaveAs("./" + 'B0_JPsiKs_DeltaTFittedOnGenMC.pdf')
665  c1.Clear()
666 
667  dotsMC = DTFrame.getHist("h_fitData")
668  dotsMCB0 = DTFrame.getHist("data_histB0")
669  dotsMCAB0 = DTFrame.getHist("data_histAB0")
670 
671  plotWithAsymmetry(
672  DTFrame,
673  self.DT,
674  dotsMC,
675  dotsMCB0,
676  dotsMCAB0,
677  curveMC,
678  curveMCB0,
679  curveMCAB0,
680  "[ps]",
681  "DeltaTFittedOnGenMC",
682  lMC,
683  textFitRes,
684  "#it{a}_{#it{CP}}(#Delta#it{t})",
685  True)
686 
687  # -----------------
688 
689  # z_tag fit
690 
691  # z_tag fit whole Delta t range
692 
693  # self.S.setVal(0.)
694  # self.A.setVal(0.)
695  self.S.setConstant(ROOT.kTRUE)
696  self.A.setConstant(ROOT.kTRUE)
697 
698  self.fitDataDeltaZtagUpsilon.Print()
699  DeltaZtagUpsilonp = ROOT.RooGenericPdf(
700  "DeltaZtagUpsilonp",
701  "DeltaZtagUpsilonp",
702  "(exp(-1*abs(DeltaZtagUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 + " +
703  "(A*(cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))) + " +
704  " S*((DM*Tau)*cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma)) - sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))))" +
705  "/(1+(Tau*DM)*(Tau*DM)))",
706  ROOT.RooArgList(
707  self.Norm,
708  self.BetaGamma,
709  self.DeltaZtagUpsilon,
710  self.Tau,
711  self.A,
712  self.S,
713  self.DM))
714  DeltaZtagUpsilonm = ROOT.RooGenericPdf(
715  "DeltaZtagUpsilonm",
716  "DeltaZtagUpsilonm",
717  "(exp(-1*abs(DeltaZtagUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 - " +
718  "(A*(cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))) + " +
719  " S*((DM*Tau)*cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma)) - sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))))" +
720  "/(1+(Tau*DM)*(Tau*DM)))",
721  ROOT.RooArgList(
722  self.Norm,
723  self.BetaGamma,
724  self.DeltaZtagUpsilon,
725  self.Tau,
726  self.A,
727  self.S,
728  self.DM))
729  DeltaZtagUpsilonModel = ROOT.RooSimultaneous("DeltaZtagUpsilonModel", "DeltaZtagUpsilonModel", self.q)
730  DeltaZtagUpsilonModel.addPdf(DeltaZtagUpsilonp, "1")
731  DeltaZtagUpsilonModel.addPdf(DeltaZtagUpsilonm, "-1")
732  DeltaZtagUpsilonModel.Print()
733 
734  fitResTag = DeltaZtagUpsilonModel.fitTo(
735  self.fitDataDeltaZtagUpsilon, ROOT.RooFit.Minos(
736  ROOT.kFALSE), ROOT.RooFit.Extended(
737  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
738  fitResTag.Print("v")
739 
740  DeltaZtagUpsilonFrame = self.DeltaZtagUpsilon.frame()
741  self.fitDataDeltaZtagUpsilon.plotOn(DeltaZtagUpsilonFrame)
742  self.fitDataDeltaZtagUpsilon.plotOn(
743  DeltaZtagUpsilonFrame,
744  ROOT.RooFit.Cut("q==q::1"),
745  ROOT.RooFit.MarkerColor(
746  ROOT.kRed),
747  ROOT.RooFit.Name("data_histB0tag"))
748  self.fitDataDeltaZtagUpsilon.plotOn(
749  DeltaZtagUpsilonFrame,
750  ROOT.RooFit.Cut("q==1"),
751  ROOT.RooFit.MarkerColor(
752  ROOT.kBlue + 1),
753  ROOT.RooFit.Name("data_histAB0tag"))
754 
755  DeltaZtagUpsilonModel.plotOn(
756  DeltaZtagUpsilonFrame,
757  ROOT.RooFit.ProjWData(
758  ROOT.RooArgSet(
759  self.q),
761  ROOT.RooFit.LineWidth(3),
762  ROOT.RooFit.LineColor(
763  ROOT.kBlack),
764  ROOT.RooFit.LineStyle(3))
765  DeltaZtagUpsilonModel.plotOn(
766  DeltaZtagUpsilonFrame,
767  ROOT.RooFit.ProjWData(
768  ROOT.RooArgSet(
769  self.q),
771  ROOT.RooFit.LineWidth(3),
772  ROOT.RooFit.LineColor(
773  ROOT.kRed),
774  ROOT.RooFit.Components("DeltaZtagUpsilonp"))
775  DeltaZtagUpsilonModel.plotOn(
776  DeltaZtagUpsilonFrame,
777  ROOT.RooFit.ProjWData(
778  ROOT.RooArgSet(
779  self.q),
781  ROOT.RooFit.LineWidth(3),
782  ROOT.RooFit.LineColor(
783  ROOT.kBlue + 1),
784  ROOT.RooFit.LineStyle(7),
785  ROOT.RooFit.Components("DeltaZtagUpsilonm"))
786  DeltaZtagUpsilonFrame.SetTitle("")
787  DeltaZtagUpsilonFrame.GetXaxis().SetTitle("(#it{z}_{tag}^{dec} - #it{z}_{tag}^{prod}){}^{gen} [mm]")
788  DeltaZtagUpsilonFrame.GetXaxis().SetTitleSize(0.05)
789  DeltaZtagUpsilonFrame.GetXaxis().SetLabelSize(0.045)
790  DeltaZtagUpsilonFrame.GetYaxis().SetTitleSize(0.05)
791  DeltaZtagUpsilonFrame.GetYaxis().SetTitleOffset(1.5)
792  DeltaZtagUpsilonFrame.GetYaxis().SetLabelSize(0.045)
793 
794  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
795  c2.cd()
796 
797  self.DeltaZtagUpsilon.Print()
798 
799  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
800  Pad.SetLeftMargin(0.15)
801  Pad.SetBottomMargin(0.15)
802  Pad.Draw()
803  Pad.cd()
804  DeltaZtagUpsilonFrame.Draw()
805 
806  curveMC = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]")
807  curveMCB0 = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]_Comp[DeltaZtagUpsilonp]")
808  curveMCAB0 = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]_Comp[DeltaZtagUpsilonm]")
809  curveMC.SetFillColor(ROOT.kWhite)
810  curveMCB0.SetFillColor(ROOT.kWhite)
811  curveMCAB0.SetFillColor(ROOT.kWhite)
812 
813  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
814  lMC.AddEntry(curveMC, "Both")
815  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
816  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
817 
818  lMC.SetTextSize(0.065)
819  lMC.Draw()
820 
821  fitResTag.Print("v")
822  c2.SaveAs("./" + 'B0_JPsiKs_DeltaZtagUpsilonFittedOnGenMC.pdf')
823  c2.Clear()
824 
825  dotsMC = DeltaZtagUpsilonFrame.getHist("h_fitDataDeltaZtagUpsilon")
826  dotsMCB0 = DeltaZtagUpsilonFrame.getHist("data_histB0tag")
827  dotsMCAB0 = DeltaZtagUpsilonFrame.getHist("data_histAB0tag")
828 
829  plotWithAsymmetry(
830  DeltaZtagUpsilonFrame,
831  self.DeltaZtagUpsilon,
832  dotsMC,
833  dotsMCB0,
834  dotsMCAB0,
835  curveMC,
836  curveMCB0,
837  curveMCAB0,
838  "[mm]",
839  "DeltaZtagUpsilonFittedOnGenMC",
840  lMC,
841  "",
842  labelAsymmetry,
843  True)
844 
845  self.S.setConstant(ROOT.kFALSE)
846  self.A.setConstant(ROOT.kFALSE)
847 
848  # -----------------
849 
850  # t_tag fit
851 
852  # t_tag fit whole Delta t range
853 
854  # self.S.setConstant(ROOT.kTRUE)
855  self.S.setVal(0.)
856  self.A.setVal(0.)
857  # self.A.setConstant(ROOT.kTRUE)
858 
859  self.fitDataTag.Print()
860  Ttagp = ROOT.RooGenericPdf(
861  "Ttagp",
862  "Ttagp",
863  "(exp(-1*abs(Ttag)/Tau)/(Norm*Tau))*(1 + " +
864  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM)))",
865  ROOT.RooArgList(
866  self.Norm,
867  self.Ttag,
868  self.Tau,
869  self.A,
870  self.S,
871  self.DM))
872  Ttagm = ROOT.RooGenericPdf(
873  "Ttagm",
874  "Ttagm",
875  "(exp(-1*abs(Ttag)/Tau)/(Norm*Tau))*(1 - " +
876  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM)))",
877  ROOT.RooArgList(
878  self.Norm,
879  self.Ttag,
880  self.Tau,
881  self.A,
882  self.S,
883  self.DM))
884  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.q)
885  TtagModel.addPdf(Ttagp, "1")
886  TtagModel.addPdf(Ttagm, "-1")
887  TtagModel.Print()
888 
889  fitResTag = TtagModel.fitTo(
890  self.fitDataTag, ROOT.RooFit.Minos(
891  ROOT.kFALSE), ROOT.RooFit.Extended(
892  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
893  fitResTag.Print("v")
894 
895  TtagFrame = self.Ttag.frame()
896  self.fitDataTag.plotOn(TtagFrame)
897  self.fitDataTag.plotOn(
898  TtagFrame,
899  ROOT.RooFit.Cut("q==q::1"),
900  ROOT.RooFit.MarkerColor(
901  ROOT.kRed),
902  ROOT.RooFit.Name("data_histB0tag"))
903  self.fitDataTag.plotOn(
904  TtagFrame,
905  ROOT.RooFit.Cut("q==1"),
906  ROOT.RooFit.MarkerColor(
907  ROOT.kBlue + 1),
908  ROOT.RooFit.Name("data_histAB0tag"))
909 
910  TtagModel.plotOn(
911  TtagFrame,
912  ROOT.RooFit.ProjWData(
913  ROOT.RooArgSet(
914  self.q),
915  self.fitDataTag),
916  ROOT.RooFit.LineWidth(3),
917  ROOT.RooFit.LineColor(
918  ROOT.kBlack),
919  ROOT.RooFit.LineStyle(3))
920  TtagModel.plotOn(
921  TtagFrame,
922  ROOT.RooFit.ProjWData(
923  ROOT.RooArgSet(
924  self.q),
925  self.fitDataTag),
926  ROOT.RooFit.LineWidth(3),
927  ROOT.RooFit.LineColor(
928  ROOT.kRed),
929  ROOT.RooFit.Components("Ttagp"))
930  TtagModel.plotOn(
931  TtagFrame,
932  ROOT.RooFit.ProjWData(
933  ROOT.RooArgSet(
934  self.q),
935  self.fitDataTag),
936  ROOT.RooFit.LineWidth(3),
937  ROOT.RooFit.LineColor(
938  ROOT.kBlue + 1),
939  ROOT.RooFit.LineStyle(7),
940  ROOT.RooFit.Components("Ttagm"))
941  TtagFrame.SetTitle("")
942  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
943  TtagFrame.GetXaxis().SetTitleSize(0.05)
944  TtagFrame.GetXaxis().SetLabelSize(0.045)
945  TtagFrame.GetYaxis().SetTitleSize(0.05)
946  TtagFrame.GetYaxis().SetTitleOffset(1.5)
947  TtagFrame.GetYaxis().SetLabelSize(0.045)
948 
949  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
950  c2.cd()
951 
952  self.Ttag.Print()
953 
954  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
955  Pad.SetLeftMargin(0.15)
956  Pad.SetBottomMargin(0.15)
957  Pad.Draw()
958  Pad.cd()
959  TtagFrame.Draw()
960 
961  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
962  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
963  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
964  curveMC.SetFillColor(ROOT.kWhite)
965  curveMCB0.SetFillColor(ROOT.kWhite)
966  curveMCAB0.SetFillColor(ROOT.kWhite)
967 
968  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
969  lMC.AddEntry(curveMC, "Both")
970  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
971  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
972 
973  lMC.SetTextSize(0.065)
974  lMC.Draw()
975 
976  fitResTag.Print("v")
977  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMC.pdf')
978  c2.Clear()
979 
980  dotsMC = TtagFrame.getHist("h_fitDataTag")
981  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
982  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
983 
984  plotWithAsymmetry(TtagFrame, self.Ttag, dotsMC, dotsMCB0, dotsMCAB0, curveMC, curveMCB0,
985  curveMCAB0, "[ps]", "TtagFittedOnGenMC", lMC, "", labelAsymmetry, True)
986 
987  # ---
988 
989  # t_tag fit positive Delta t range
990 
991  # self.S.setConstant(ROOT.kTRUE)
992  # self.S.setVal(0.)
993  # self.A.setVal(0.)
994  # self.A.setConstant(ROOT.kTRUE)
995 
996  self.fitDataTagDeltaTPos.Print()
997  Ttagp = ROOT.RooGenericPdf(
998  "Ttagp",
999  "Ttagp",
1000  "(exp(-2*abs(Ttag)/Tau)/(Tau))*(1 + (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1001  ROOT.RooArgList(
1002  self.Norm,
1003  self.Ttag,
1004  self.Tau,
1005  self.A,
1006  self.S,
1007  self.DM))
1008  Ttagm = ROOT.RooGenericPdf(
1009  "Ttagm",
1010  "Ttagm",
1011  "(exp(-2*abs(Ttag)/Tau)/(Tau))*(1 - (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1012  ROOT.RooArgList(
1013  self.Norm,
1014  self.Ttag,
1015  self.Tau,
1016  self.A,
1017  self.S,
1018  self.DM))
1019  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.q)
1020  TtagModel.addPdf(Ttagp, "1")
1021  TtagModel.addPdf(Ttagm, "-1")
1022  TtagModel.Print()
1023 
1024  fitResTagPos = TtagModel.fitTo(
1025  self.fitDataTagDeltaTPos, ROOT.RooFit.Minos(
1026  ROOT.kFALSE), ROOT.RooFit.Extended(
1027  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1028  fitResTagPos.Print("v")
1029 
1030  TtagFrame = self.TtagPosNeg.frame()
1031  self.fitDataTagDeltaTPos.plotOn(TtagFrame, ROOT.RooFit.Cut("Ttag < " + str(TtagFrame.GetXaxis().GetXmax())))
1032  self.fitDataTagDeltaTPos.plotOn(TtagFrame,
1033  ROOT.RooFit.Cut("q==q::1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1034  ROOT.RooFit.MarkerColor(ROOT.kRed),
1035  ROOT.RooFit.Name("data_histB0tag"))
1036  self.fitDataTagDeltaTPos.plotOn(TtagFrame,
1037  ROOT.RooFit.Cut("q==1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1038  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1039  ROOT.RooFit.Name("data_histAB0tag"))
1040 
1041  TtagModel.plotOn(TtagFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q), self.fitDataTagDeltaTPos),
1042  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1043  TtagModel.plotOn(
1044  TtagFrame,
1045  ROOT.RooFit.ProjWData(
1046  ROOT.RooArgSet(
1047  self.q),
1048  self.fitDataTagDeltaTPos),
1049  ROOT.RooFit.LineWidth(3),
1050  ROOT.RooFit.LineColor(
1051  ROOT.kRed),
1052  ROOT.RooFit.Components("Ttagp"))
1053  TtagModel.plotOn(
1054  TtagFrame,
1055  ROOT.RooFit.ProjWData(
1056  ROOT.RooArgSet(
1057  self.q),
1058  self.fitDataTagDeltaTPos),
1059  ROOT.RooFit.LineWidth(3),
1060  ROOT.RooFit.LineColor(
1061  ROOT.kBlue + 1),
1062  ROOT.RooFit.LineStyle(7),
1063  ROOT.RooFit.Components("Ttagm"))
1064  TtagFrame.SetTitle("")
1065  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
1066  TtagFrame.GetXaxis().SetTitleSize(0.05)
1067  TtagFrame.GetXaxis().SetLabelSize(0.045)
1068  TtagFrame.GetYaxis().SetTitleSize(0.05)
1069  TtagFrame.GetYaxis().SetTitleOffset(1.5)
1070  TtagFrame.GetYaxis().SetLabelSize(0.045)
1071 
1072  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1073  c2.cd()
1074 
1075  self.Ttag.Print()
1076 
1077  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1078  Pad.SetLeftMargin(0.15)
1079  Pad.SetBottomMargin(0.15)
1080  Pad.Draw()
1081  Pad.cd()
1082  TtagFrame.Draw()
1083 
1084  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
1085  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
1086  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
1087  curveMC.SetFillColor(ROOT.kWhite)
1088  curveMCB0.SetFillColor(ROOT.kWhite)
1089  curveMCAB0.SetFillColor(ROOT.kWhite)
1090 
1091  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1092  lMC.AddEntry(curveMC, "Both")
1093  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1094  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1095 
1096  lMC.SetTextSize(0.065)
1097  lMC.Draw()
1098 
1099  fitResTagPos.Print("v")
1100  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMCPositiveDeltaT.pdf')
1101  c2.Clear()
1102 
1103  dotsMC = TtagFrame.getHist("h_fitDataTag")
1104  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
1105  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
1106 
1107  plotWithAsymmetry(
1108  TtagFrame,
1109  self.TtagPosNeg,
1110  dotsMC,
1111  dotsMCB0,
1112  dotsMCAB0,
1113  curveMC,
1114  curveMCB0,
1115  curveMCAB0,
1116  "[ps]",
1117  "TtagFittedOnGenMCPositiveDeltaT",
1118  lMC,
1119  "",
1120  labelAsymmetry,
1121  True)
1122 
1123  # t_tag fit negative Delta t range
1124 
1125  # self.S.setConstant(ROOT.kTRUE)
1126  self.S.setVal(0.)
1127  self.A.setVal(0.)
1128  # self.A.setConstant(ROOT.kTRUE)
1129  self.A.setConstant(ROOT.kFALSE)
1130 
1131  self.fitDataTagDeltaTNeg.Print()
1132  Ttagp = ROOT.RooGenericPdf(
1133  "Ttagp",
1134  "Ttagp",
1135  "(exp(-1*abs(Ttag)/Tau)*(1 + " +
1136  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM))) - " +
1137  "(exp(-2*abs(Ttag)/Tau))*(1 + (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1138  ROOT.RooArgList(
1139  self.Ttag,
1140  self.Tau,
1141  self.A,
1142  self.S,
1143  self.DM))
1144  Ttagm = ROOT.RooGenericPdf(
1145  "Ttagm",
1146  "Ttagm",
1147  "(exp(-1*abs(Ttag)/Tau)*(1 - " +
1148  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM))) - " +
1149  "(exp(-2*abs(Ttag)/Tau))*(1 - (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1150  ROOT.RooArgList(
1151  self.Ttag,
1152  self.Tau,
1153  self.A,
1154  self.S,
1155  self.DM))
1156  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.q)
1157  TtagModel.addPdf(Ttagp, "1")
1158  TtagModel.addPdf(Ttagm, "-1")
1159  TtagModel.Print()
1160 
1161  fitResTagNeg = TtagModel.fitTo(
1162  self.fitDataTagDeltaTNeg, ROOT.RooFit.Minos(
1163  ROOT.kFALSE), ROOT.RooFit.Extended(
1164  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1165  fitResTagNeg.Print("v")
1166 
1167  TtagFrame = self.TtagPosNeg.frame()
1168  self.fitDataTagDeltaTNeg.plotOn(TtagFrame, ROOT.RooFit.Cut("Ttag < " + str(TtagFrame.GetXaxis().GetXmax())))
1169  self.fitDataTagDeltaTNeg.plotOn(TtagFrame,
1170  ROOT.RooFit.Cut("q==q::1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1171  ROOT.RooFit.MarkerColor(ROOT.kRed),
1172  ROOT.RooFit.Name("data_histB0tag"))
1173  self.fitDataTagDeltaTNeg.plotOn(TtagFrame,
1174  ROOT.RooFit.Cut("q==1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1175  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1176  ROOT.RooFit.Name("data_histAB0tag"))
1177 
1178  TtagModel.plotOn(TtagFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q), self.fitDataTagDeltaTNeg),
1179  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1180  TtagModel.plotOn(
1181  TtagFrame,
1182  ROOT.RooFit.ProjWData(
1183  ROOT.RooArgSet(
1184  self.q),
1185  self.fitDataTagDeltaTNeg),
1186  ROOT.RooFit.LineWidth(3),
1187  ROOT.RooFit.LineColor(
1188  ROOT.kRed),
1189  ROOT.RooFit.Components("Ttagp"))
1190  TtagModel.plotOn(
1191  TtagFrame,
1192  ROOT.RooFit.ProjWData(
1193  ROOT.RooArgSet(
1194  self.q),
1195  self.fitDataTagDeltaTNeg),
1196  ROOT.RooFit.LineWidth(3),
1197  ROOT.RooFit.LineColor(
1198  ROOT.kBlue + 1),
1199  ROOT.RooFit.LineStyle(7),
1200  ROOT.RooFit.Components("Ttagm"))
1201  TtagFrame.SetTitle("")
1202  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
1203  TtagFrame.GetXaxis().SetTitleSize(0.05)
1204  TtagFrame.GetXaxis().SetLabelSize(0.045)
1205  TtagFrame.GetYaxis().SetTitleSize(0.05)
1206  TtagFrame.GetYaxis().SetTitleOffset(1.5)
1207  TtagFrame.GetYaxis().SetLabelSize(0.045)
1208 
1209  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1210  c2.cd()
1211 
1212  self.Ttag.Print()
1213 
1214  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1215  Pad.SetLeftMargin(0.15)
1216  Pad.SetBottomMargin(0.15)
1217  Pad.Draw()
1218  Pad.cd()
1219  TtagFrame.Draw()
1220  TtagFrame.Print()
1221 
1222  # lMC = ROOT.TLegend(0.7, 0.63, 0.9, 0.89)
1223 
1224  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
1225  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
1226  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
1227  curveMC.SetFillColor(ROOT.kWhite)
1228  curveMCB0.SetFillColor(ROOT.kWhite)
1229  curveMCAB0.SetFillColor(ROOT.kWhite)
1230 
1231  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1232  lMC.AddEntry(curveMC, "Both")
1233  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1234  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1235 
1236  lMC.SetTextSize(0.065)
1237  lMC.Draw()
1238 
1239  fitResTagNeg.Print("v")
1240  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMCNegativeDeltaT.pdf')
1241  c2.Clear()
1242 
1243  dotsMC = TtagFrame.getHist("h_fitDataTag")
1244  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
1245  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
1246 
1247  plotWithAsymmetry(
1248  TtagFrame,
1249  self.TtagPosNeg,
1250  dotsMC,
1251  dotsMCB0,
1252  dotsMCAB0,
1253  curveMC,
1254  curveMCB0,
1255  curveMCAB0,
1256  "[ps]",
1257  "TtagFittedOnGenMCNegativeDeltaT",
1258  lMC,
1259  "",
1260  labelAsymmetry,
1261  True)
1262 
1263  # ----------------------
1264 
1265  # z_sig fit
1266 
1267  # z_sig fit whole Delta t range
1268 
1269  # self.S.setVal(0.)
1270  # self.A.setVal(0.)
1271  self.S.setConstant(ROOT.kTRUE)
1272  self.A.setConstant(ROOT.kTRUE)
1273 
1274  self.fitDataDeltaZsigUpsilon.Print()
1275  DeltaZsigUpsilonp = ROOT.RooGenericPdf(
1276  "DeltaZsigUpsilonp",
1277  "DeltaZsigUpsilonp",
1278  "(exp(-1*abs(DeltaZsigUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 + " +
1279  "(A*(cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))) + " +
1280  "S*(-(DM*Tau)*cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma)) + sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))))" +
1281  "/(1+(Tau*DM)*(Tau*DM)))",
1282  ROOT.RooArgList(
1283  self.Norm,
1284  self.BetaGamma,
1285  self.DeltaZsigUpsilon,
1286  self.Tau,
1287  self.A,
1288  self.S,
1289  self.DM))
1290  DeltaZsigUpsilonm = ROOT.RooGenericPdf(
1291  "DeltaZsigUpsilonm",
1292  "DeltaZsigUpsilonm",
1293  "(exp(-1*abs(DeltaZsigUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 - "
1294  "(A*(cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))) + "
1295  "S*(-(DM*Tau)*cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma)) + sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))))"
1296  "/(1+(Tau*DM)*(Tau*DM)))",
1297  ROOT.RooArgList(
1298  self.Norm,
1299  self.BetaGamma,
1300  self.DeltaZsigUpsilon,
1301  self.Tau,
1302  self.A,
1303  self.S,
1304  self.DM))
1305  DeltaZsigUpsilonModel = ROOT.RooSimultaneous("DeltaZsigUpsilonModel", "DeltaZsigUpsilonModel", self.q)
1306  DeltaZsigUpsilonModel.addPdf(DeltaZsigUpsilonp, "1")
1307  DeltaZsigUpsilonModel.addPdf(DeltaZsigUpsilonm, "-1")
1308  DeltaZsigUpsilonModel.Print()
1309 
1310  fitResTag = DeltaZsigUpsilonModel.fitTo(
1311  self.fitDataDeltaZsigUpsilon, ROOT.RooFit.Minos(
1312  ROOT.kFALSE), ROOT.RooFit.Extended(
1313  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1314  fitResTag.Print("v")
1315 
1316  DeltaZsigUpsilonFrame = self.DeltaZsigUpsilon.frame()
1317  self.fitDataDeltaZsigUpsilon.plotOn(DeltaZsigUpsilonFrame)
1318  self.fitDataDeltaZsigUpsilon.plotOn(
1319  DeltaZsigUpsilonFrame,
1320  ROOT.RooFit.Cut("q==q::1"),
1321  ROOT.RooFit.MarkerColor(
1322  ROOT.kRed),
1323  ROOT.RooFit.Name("data_histB0sig"))
1324  self.fitDataDeltaZsigUpsilon.plotOn(
1325  DeltaZsigUpsilonFrame,
1326  ROOT.RooFit.Cut("q==1"),
1327  ROOT.RooFit.MarkerColor(
1328  ROOT.kBlue + 1),
1329  ROOT.RooFit.Name("data_histAB0sig"))
1330 
1331  DeltaZsigUpsilonModel.plotOn(
1332  DeltaZsigUpsilonFrame,
1333  ROOT.RooFit.ProjWData(
1334  ROOT.RooArgSet(
1335  self.q),
1337  ROOT.RooFit.LineWidth(3),
1338  ROOT.RooFit.LineColor(
1339  ROOT.kBlack),
1340  ROOT.RooFit.LineStyle(3))
1341  DeltaZsigUpsilonModel.plotOn(
1342  DeltaZsigUpsilonFrame,
1343  ROOT.RooFit.ProjWData(
1344  ROOT.RooArgSet(
1345  self.q),
1347  ROOT.RooFit.LineWidth(3),
1348  ROOT.RooFit.LineColor(
1349  ROOT.kRed),
1350  ROOT.RooFit.Components("DeltaZsigUpsilonp"))
1351  DeltaZsigUpsilonModel.plotOn(
1352  DeltaZsigUpsilonFrame,
1353  ROOT.RooFit.ProjWData(
1354  ROOT.RooArgSet(
1355  self.q),
1357  ROOT.RooFit.LineWidth(3),
1358  ROOT.RooFit.LineColor(
1359  ROOT.kBlue + 1),
1360  ROOT.RooFit.LineStyle(7),
1361  ROOT.RooFit.Components("DeltaZsigUpsilonm"))
1362  DeltaZsigUpsilonFrame.SetTitle("")
1363  DeltaZsigUpsilonFrame.GetXaxis().SetTitle("(#it{z}_{sig}^{dec} - #it{z}_{sig}^{prod}){}^{gen} [mm]")
1364  DeltaZsigUpsilonFrame.GetXaxis().SetTitleSize(0.05)
1365  DeltaZsigUpsilonFrame.GetXaxis().SetLabelSize(0.045)
1366  DeltaZsigUpsilonFrame.GetYaxis().SetTitleSize(0.05)
1367  DeltaZsigUpsilonFrame.GetYaxis().SetTitleOffset(1.5)
1368  DeltaZsigUpsilonFrame.GetYaxis().SetLabelSize(0.045)
1369 
1370  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1371  c2.cd()
1372 
1373  self.DeltaZsigUpsilon.Print()
1374 
1375  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1376  Pad.SetLeftMargin(0.15)
1377  Pad.SetBottomMargin(0.15)
1378  Pad.Draw()
1379  Pad.cd()
1380  DeltaZsigUpsilonFrame.Draw()
1381 
1382  curveMC = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]")
1383  curveMCB0 = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]_Comp[DeltaZsigUpsilonp]")
1384  curveMCAB0 = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]_Comp[DeltaZsigUpsilonm]")
1385  curveMC.SetFillColor(ROOT.kWhite)
1386  curveMCB0.SetFillColor(ROOT.kWhite)
1387  curveMCAB0.SetFillColor(ROOT.kWhite)
1388 
1389  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1390  lMC.AddEntry(curveMC, "Both")
1391  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1392  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1393 
1394  lMC.SetTextSize(0.065)
1395  lMC.Draw()
1396 
1397  fitResTag.Print("v")
1398  c2.SaveAs("./" + 'B0_JPsiKs_DeltaZsigUpsilonFittedOnGenMC.pdf')
1399  c2.Clear()
1400 
1401  dotsMC = DeltaZsigUpsilonFrame.getHist("h_fitDataDeltaZsigUpsilon")
1402  dotsMCB0 = DeltaZsigUpsilonFrame.getHist("data_histB0sig")
1403  dotsMCAB0 = DeltaZsigUpsilonFrame.getHist("data_histAB0sig")
1404 
1405  plotWithAsymmetry(
1406  DeltaZsigUpsilonFrame,
1407  self.DeltaZsigUpsilon,
1408  dotsMC,
1409  dotsMCB0,
1410  dotsMCAB0,
1411  curveMC,
1412  curveMCB0,
1413  curveMCAB0,
1414  "[mm]",
1415  "DeltaZsigUpsilonFittedOnGenMC",
1416  lMC,
1417  "",
1418  labelAsymmetry,
1419  True)
1420 
1421  self.S.setConstant(ROOT.kFALSE)
1422  self.A.setConstant(ROOT.kFALSE)
1423 
1424  # --- Signal time fit --------
1425 
1426  self.fitDataSig.Print()
1427  Tsigp = ROOT.RooGenericPdf(
1428  "Tsigp",
1429  "Tsigp",
1430  "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))*(1 + " +
1431  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM)))",
1432  ROOT.RooArgList(
1433  self.Norm,
1434  self.Tsig,
1435  self.Tau,
1436  self.A,
1437  self.S,
1438  self.DM))
1439  Tsigm = ROOT.RooGenericPdf(
1440  "Tsigm",
1441  "Tsigm",
1442  "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))*(1 - " +
1443  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM)))",
1444  ROOT.RooArgList(
1445  self.Norm,
1446  self.Tsig,
1447  self.Tau,
1448  self.A,
1449  self.S,
1450  self.DM))
1451  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.q)
1452  TsigModel.addPdf(Tsigp, "1")
1453  TsigModel.addPdf(Tsigm, "-1")
1454  TsigModel.Print()
1455 
1456  fitResSig = TsigModel.fitTo(
1457  self.fitDataSig, ROOT.RooFit.Minos(
1458  ROOT.kFALSE), ROOT.RooFit.Extended(
1459  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1460  fitResSig.Print("v")
1461 
1462  TsigFrame = self.Tsig.frame()
1463  self.fitDataSig.plotOn(TsigFrame)
1464  self.fitDataSig.plotOn(
1465  TsigFrame,
1466  ROOT.RooFit.Cut("q==q::1"),
1467  ROOT.RooFit.MarkerColor(
1468  ROOT.kRed),
1469  ROOT.RooFit.Name("data_histB0sig"))
1470  self.fitDataSig.plotOn(
1471  TsigFrame,
1472  ROOT.RooFit.Cut("q==1"),
1473  ROOT.RooFit.MarkerColor(
1474  ROOT.kBlue + 1),
1475  ROOT.RooFit.Name("data_histAB0sig"))
1476 
1477  TsigModel.plotOn(
1478  TsigFrame,
1479  ROOT.RooFit.ProjWData(
1480  ROOT.RooArgSet(
1481  self.qsig),
1482  self.fitDataSig),
1483  ROOT.RooFit.LineWidth(3),
1484  ROOT.RooFit.LineColor(
1485  ROOT.kBlack),
1486  ROOT.RooFit.LineStyle(3))
1487  TsigModel.plotOn(
1488  TsigFrame,
1489  ROOT.RooFit.ProjWData(
1490  ROOT.RooArgSet(
1491  self.qsig),
1492  self.fitDataSig),
1493  ROOT.RooFit.LineWidth(3),
1494  ROOT.RooFit.LineColor(
1495  ROOT.kRed),
1496  ROOT.RooFit.Components("Tsigp"))
1497  TsigModel.plotOn(
1498  TsigFrame,
1499  ROOT.RooFit.ProjWData(
1500  ROOT.RooArgSet(
1501  self.qsig),
1502  self.fitDataSig),
1503  ROOT.RooFit.LineWidth(3),
1504  ROOT.RooFit.LineColor(
1505  ROOT.kBlue + 1),
1506  ROOT.RooFit.LineStyle(7),
1507  ROOT.RooFit.Components("Tsigm"))
1508  TsigFrame.SetTitle("")
1509  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1510  TsigFrame.GetXaxis().SetTitleSize(0.05)
1511  TsigFrame.GetXaxis().SetLabelSize(0.045)
1512  TsigFrame.GetYaxis().SetTitleSize(0.05)
1513  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1514  TsigFrame.GetYaxis().SetLabelSize(0.045)
1515 
1516  c3 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1517  c3.cd()
1518 
1519  self.Tsig.Print()
1520 
1521  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1522  Pad.SetLeftMargin(0.15)
1523  Pad.SetBottomMargin(0.15)
1524  Pad.Draw()
1525  Pad.cd()
1526  TsigFrame.Draw()
1527 
1528  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1529  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1530  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1531  curveMC.SetFillColor(ROOT.kWhite)
1532  curveMCB0.SetFillColor(ROOT.kWhite)
1533  curveMCAB0.SetFillColor(ROOT.kWhite)
1534 
1535  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1536  lMC.AddEntry(curveMC, "Both")
1537  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1538  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1539 
1540  lMC.SetTextSize(0.065)
1541  lMC.Draw()
1542 
1543  fitResSig.Print("v")
1544  c3.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMC.pdf')
1545  c3.Clear()
1546 
1547  dotsMC = TsigFrame.getHist("h_fitDataSig")
1548  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1549  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1550 
1551  plotWithAsymmetry(TsigFrame, self.Tsig, dotsMC, dotsMCB0, dotsMCAB0, curveMC, curveMCB0,
1552  curveMCAB0, "[ps]", "TsigFittedOnGenMC", lMC, "", labelAsymmetry, True)
1553 
1554  # --- ------
1555 
1556  # t_sig fit positive Delta t range
1557 
1558  # self.S.setConstant(ROOT.kTRUE)
1559  # self.S.setVal(0.)
1560  # self.A.setVal(0.)
1561  # self.A.setConstant(ROOT.kTRUE)
1562 
1563  self.fitDataSigDeltaTPos.Print()
1564 
1565  Tsigp = ROOT.RooGenericPdf(
1566  "Tsigp",
1567  "Tsigp",
1568  "(exp(-1*abs(Tsig)/Tau)*(1 + " +
1569  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM))) - " +
1570  "(exp(-2*abs(Tsig)/Tau))*(1 + (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1571  ROOT.RooArgList(
1572  self.Tsig,
1573  self.Tau,
1574  self.A,
1575  self.S,
1576  self.DM))
1577  Tsigm = ROOT.RooGenericPdf(
1578  "Tsigm",
1579  "Tsigm",
1580  "(exp(-1*abs(Tsig)/Tau)*(1 - " +
1581  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM))) - " +
1582  "(exp(-2*abs(Tsig)/Tau))*(1 - (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1583  ROOT.RooArgList(
1584  self.Tsig,
1585  self.Tau,
1586  self.A,
1587  self.S,
1588  self.DM))
1589 
1590  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.q)
1591  TsigModel.addPdf(Tsigp, "1")
1592  TsigModel.addPdf(Tsigm, "-1")
1593  TsigModel.Print()
1594 
1595  fitResSigPos = TsigModel.fitTo(
1596  self.fitDataSigDeltaTPos, ROOT.RooFit.Minos(
1597  ROOT.kFALSE), ROOT.RooFit.Extended(
1598  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1599  fitResSigPos.Print("v")
1600 
1601  TsigFrame = self.TsigPosNeg.frame()
1602  self.fitDataSigDeltaTPos.plotOn(TsigFrame, ROOT.RooFit.Cut("Tsig < " + str(TsigFrame.GetXaxis().GetXmax())))
1603  self.fitDataSigDeltaTPos.plotOn(TsigFrame,
1604  ROOT.RooFit.Cut("q==q::1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1605  ROOT.RooFit.MarkerColor(ROOT.kRed),
1606  ROOT.RooFit.Name("data_histB0sig"))
1607  self.fitDataSigDeltaTPos.plotOn(TsigFrame,
1608  ROOT.RooFit.Cut("q==1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1609  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1610  ROOT.RooFit.Name("data_histAB0sig"))
1611 
1612  TsigModel.plotOn(TsigFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q), self.fitDataSigDeltaTPos),
1613  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1614  TsigModel.plotOn(
1615  TsigFrame,
1616  ROOT.RooFit.ProjWData(
1617  ROOT.RooArgSet(
1618  self.q),
1619  self.fitDataSigDeltaTPos),
1620  ROOT.RooFit.LineWidth(3),
1621  ROOT.RooFit.LineColor(
1622  ROOT.kRed),
1623  ROOT.RooFit.Components("Tsigp"))
1624  TsigModel.plotOn(
1625  TsigFrame,
1626  ROOT.RooFit.ProjWData(
1627  ROOT.RooArgSet(
1628  self.q),
1629  self.fitDataSigDeltaTPos),
1630  ROOT.RooFit.LineWidth(3),
1631  ROOT.RooFit.LineColor(
1632  ROOT.kBlue + 1),
1633  ROOT.RooFit.LineStyle(7),
1634  ROOT.RooFit.Components("Tsigm"))
1635  TsigFrame.SetTitle("")
1636  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1637  TsigFrame.GetXaxis().SetTitleSize(0.05)
1638  TsigFrame.GetXaxis().SetLabelSize(0.045)
1639  TsigFrame.GetYaxis().SetTitleSize(0.05)
1640  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1641  TsigFrame.GetYaxis().SetLabelSize(0.045)
1642 
1643  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1644  c2.cd()
1645 
1646  self.Tsig.Print()
1647 
1648  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1649  Pad.SetLeftMargin(0.15)
1650  Pad.SetBottomMargin(0.15)
1651  Pad.Draw()
1652  Pad.cd()
1653  TsigFrame.Draw()
1654 
1655  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1656  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1657  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1658  curveMC.SetFillColor(ROOT.kWhite)
1659  curveMCB0.SetFillColor(ROOT.kWhite)
1660  curveMCAB0.SetFillColor(ROOT.kWhite)
1661 
1662  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1663  lMC.AddEntry(curveMC, "Both")
1664  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1665  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1666 
1667  lMC.SetTextSize(0.065)
1668  lMC.Draw()
1669 
1670  fitResSigPos.Print("v")
1671  c2.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMCPositiveDeltaT.pdf')
1672  c2.Clear()
1673 
1674  dotsMC = TsigFrame.getHist("h_fitDataSig")
1675  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1676  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1677 
1678  plotWithAsymmetry(
1679  TsigFrame,
1680  self.TsigPosNeg,
1681  dotsMC,
1682  dotsMCB0,
1683  dotsMCAB0,
1684  curveMC,
1685  curveMCB0,
1686  curveMCAB0,
1687  "[ps]",
1688  "TsigFittedOnGenMCPositiveDeltaT",
1689  lMC,
1690  "",
1691  labelAsymmetry,
1692  True)
1693 
1694  # t_sig fit negative Delta t range
1695 
1696  # self.S.setConstant(ROOT.kTRUE)
1697  self.S.setVal(0.)
1698  self.A.setVal(0.)
1699  # self.A.setConstant(ROOT.kTRUE)
1700  self.A.setConstant(ROOT.kFALSE)
1701 
1702  self.fitDataSigDeltaTNeg.Print()
1703 
1704  Tsigp = ROOT.RooGenericPdf(
1705  "Tsigp",
1706  "Tsigp",
1707  "(exp(-2*abs(Tsig)/Tau)/(Tau))*(1 + (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1708  ROOT.RooArgList(
1709  self.Norm,
1710  self.Tsig,
1711  self.Tau,
1712  self.A,
1713  self.S,
1714  self.DM))
1715  Tsigm = ROOT.RooGenericPdf(
1716  "Tsigm",
1717  "Tsigm",
1718  "(exp(-2*abs(Tsig)/Tau)/(Tau))*(1 - (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1719  ROOT.RooArgList(
1720  self.Norm,
1721  self.Tsig,
1722  self.Tau,
1723  self.A,
1724  self.S,
1725  self.DM))
1726 
1727  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.q)
1728  TsigModel.addPdf(Tsigp, "1")
1729  TsigModel.addPdf(Tsigm, "-1")
1730  TsigModel.Print()
1731 
1732  fitResSigNeg = TsigModel.fitTo(
1733  self.fitDataSigDeltaTNeg, ROOT.RooFit.Minos(
1734  ROOT.kFALSE), ROOT.RooFit.Extended(
1735  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1736  fitResSigNeg.Print("v")
1737 
1738  TsigFrame = self.TsigPosNeg.frame()
1739  self.fitDataSigDeltaTNeg.plotOn(TsigFrame, ROOT.RooFit.Cut("Tsig < " + str(TsigFrame.GetXaxis().GetXmax())))
1740  self.fitDataSigDeltaTNeg.plotOn(TsigFrame,
1741  ROOT.RooFit.Cut("q==q::1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1742  ROOT.RooFit.MarkerColor(ROOT.kRed),
1743  ROOT.RooFit.Name("data_histB0sig"))
1744  self.fitDataSigDeltaTNeg.plotOn(TsigFrame,
1745  ROOT.RooFit.Cut("q==1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1746  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1747  ROOT.RooFit.Name("data_histAB0sig"))
1748 
1749  TsigModel.plotOn(TsigFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q), self.fitDataSigDeltaTNeg),
1750  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1751  TsigModel.plotOn(
1752  TsigFrame,
1753  ROOT.RooFit.ProjWData(
1754  ROOT.RooArgSet(
1755  self.q),
1756  self.fitDataSigDeltaTNeg),
1757  ROOT.RooFit.LineWidth(3),
1758  ROOT.RooFit.LineColor(
1759  ROOT.kRed),
1760  ROOT.RooFit.Components("Tsigp"))
1761  TsigModel.plotOn(
1762  TsigFrame,
1763  ROOT.RooFit.ProjWData(
1764  ROOT.RooArgSet(
1765  self.q),
1766  self.fitDataSigDeltaTNeg),
1767  ROOT.RooFit.LineWidth(3),
1768  ROOT.RooFit.LineColor(
1769  ROOT.kBlue + 1),
1770  ROOT.RooFit.LineStyle(7),
1771  ROOT.RooFit.Components("Tsigm"))
1772  TsigFrame.SetTitle("")
1773  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1774  TsigFrame.GetXaxis().SetTitleSize(0.05)
1775  TsigFrame.GetXaxis().SetLabelSize(0.045)
1776  TsigFrame.GetYaxis().SetTitleSize(0.05)
1777  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1778  TsigFrame.GetYaxis().SetLabelSize(0.045)
1779 
1780  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1781  c2.cd()
1782 
1783  self.Tsig.Print()
1784 
1785  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1786  Pad.SetLeftMargin(0.15)
1787  Pad.SetBottomMargin(0.15)
1788  Pad.Draw()
1789  Pad.cd()
1790  TsigFrame.Draw()
1791  TsigFrame.Print()
1792 
1793  # lMC = ROOT.TLegend(0.7, 0.63, 0.9, 0.89)
1794 
1795  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1796  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1797  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1798  curveMC.SetFillColor(ROOT.kWhite)
1799  curveMCB0.SetFillColor(ROOT.kWhite)
1800  curveMCAB0.SetFillColor(ROOT.kWhite)
1801 
1802  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1803  lMC.AddEntry(curveMC, "Both")
1804  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1805  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1806 
1807  lMC.SetTextSize(0.065)
1808  lMC.Draw()
1809 
1810  fitResSigNeg.Print("v")
1811  c2.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMCNegativeDeltaT.pdf')
1812  c2.Clear()
1813 
1814  dotsMC = TsigFrame.getHist("h_fitDataSig")
1815  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1816  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1817 
1818  plotWithAsymmetry(
1819  TsigFrame,
1820  self.TsigPosNeg,
1821  dotsMC,
1822  dotsMCB0,
1823  dotsMCAB0,
1824  curveMC,
1825  curveMCB0,
1826  curveMCAB0,
1827  "[ps]",
1828  "TsigFittedOnGenMCNegativeDeltaT",
1829  lMC,
1830  "",
1831  labelAsymmetry,
1832  True)
1833 
1834  # ----- Signal time fit qsig -------
1835 
1836  self.Tau.setConstant(ROOT.kFALSE)
1837 
1838  self.fitDataSig.Print()
1839  Tsigp = ROOT.RooGenericPdf("Tsigp", "Tsigp", "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))",
1840  ROOT.RooArgList(self.Norm, self.Tsig, self.Tau, self.DM))
1841  Tsigm = ROOT.RooGenericPdf("Tsigm", "Tsigm", "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))",
1842  ROOT.RooArgList(self.Norm, self.Tsig, self.Tau, self.DM))
1843  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.qsig)
1844  TsigModel.addPdf(Tsigp, "1")
1845  TsigModel.addPdf(Tsigm, "-1")
1846  TsigModel.Print()
1847 
1848  fitResSig = TsigModel.fitTo(
1849  self.fitDataSig, ROOT.RooFit.Minos(
1850  ROOT.kFALSE), ROOT.RooFit.Extended(
1851  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1852  fitResSig.Print("v")
1853 
1854  TsigFrame = self.Tsig.frame()
1855  self.fitDataSig.plotOn(TsigFrame)
1856  self.fitDataSig.plotOn(
1857  TsigFrame,
1858  ROOT.RooFit.Cut("qsig==qsig::1"),
1859  ROOT.RooFit.MarkerColor(
1860  ROOT.kRed),
1861  ROOT.RooFit.Name("data_histB0sig"))
1862  self.fitDataSig.plotOn(
1863  TsigFrame,
1864  ROOT.RooFit.Cut("qsig==1"),
1865  ROOT.RooFit.MarkerColor(
1866  ROOT.kBlue + 1),
1867  ROOT.RooFit.Name("data_histAB0sig"))
1868 
1869  TsigModel.plotOn(
1870  TsigFrame,
1871  ROOT.RooFit.ProjWData(
1872  ROOT.RooArgSet(
1873  self.qsig),
1874  self.fitDataSig),
1875  ROOT.RooFit.LineWidth(3),
1876  ROOT.RooFit.LineColor(
1877  ROOT.kBlack),
1878  ROOT.RooFit.LineStyle(3))
1879  TsigModel.plotOn(
1880  TsigFrame,
1881  ROOT.RooFit.ProjWData(
1882  ROOT.RooArgSet(
1883  self.qsig),
1884  self.fitDataSig),
1885  ROOT.RooFit.LineWidth(3),
1886  ROOT.RooFit.LineColor(
1887  ROOT.kRed),
1888  ROOT.RooFit.Components("Tsigp"))
1889  TsigModel.plotOn(
1890  TsigFrame,
1891  ROOT.RooFit.ProjWData(
1892  ROOT.RooArgSet(
1893  self.qsig),
1894  self.fitDataSig),
1895  ROOT.RooFit.LineWidth(3),
1896  ROOT.RooFit.LineColor(
1897  ROOT.kBlue + 1),
1898  ROOT.RooFit.LineStyle(7),
1899  ROOT.RooFit.Components("Tsigm"))
1900  TsigFrame.SetTitle("")
1901  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1902  TsigFrame.GetXaxis().SetTitleSize(0.05)
1903  TsigFrame.GetXaxis().SetLabelSize(0.045)
1904  TsigFrame.GetYaxis().SetTitleSize(0.05)
1905  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1906  TsigFrame.GetYaxis().SetLabelSize(0.045)
1907 
1908  c3 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1909  c3.cd()
1910 
1911  self.Tsig.Print()
1912 
1913  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1914  Pad.SetLeftMargin(0.15)
1915  Pad.SetBottomMargin(0.15)
1916  Pad.Draw()
1917  Pad.cd()
1918  TsigFrame.Draw()
1919 
1920  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1921  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1922  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1923  curveMC.SetFillColor(ROOT.kWhite)
1924  curveMCB0.SetFillColor(ROOT.kWhite)
1925  curveMCAB0.SetFillColor(ROOT.kWhite)
1926 
1927  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1928  lMC.AddEntry(curveMC, "Both")
1929  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1930  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1931 
1932  lMC.SetTextSize(0.065)
1933  lMC.Draw()
1934 
1935  fitResSig.Print("v")
1936  c3.SaveAs("./" + 'B0_JPsiKs_TsigQsigFittedOnGenMC.pdf')
1937  c3.Clear()
1938 
1939  B0tagPercentage = self.B0sInTagSide / (self.B0barsInTagSide + self.B0sInTagSide) * 100
1940  B0tagbarsPercentage = self.B0barsInTagSide / (self.B0barsInTagSide + self.B0sInTagSide) * 100
1941  Diff = B0tagbarsPercentage - B0tagPercentage
1942  print("B0 Tag Flavor: ")
1943  print(
1944  "B0s = ",
1945  " %0.4f" %
1946  B0tagPercentage,
1947  "% B0bars = ",
1948  " %0.4f" %
1949  B0tagbarsPercentage,
1950  "% Diff = ",
1951  " %0.2f" %
1952  Diff,
1953  "%")
1954  B0signalPercentage = self.B0sInSignalSide / (self.B0barsInSignalSide + self.B0sInSignalSide) * 100
1955  B0signalbarsPercentage = self.B0barsInSignalSide / (self.B0barsInSignalSide + self.B0sInSignalSide) * 100
1956  Diff = B0signalbarsPercentage - B0signalPercentage
1957  print("B0 Signal Flavor: ")
1958  print(
1959  "B0s = ",
1960  " %0.4f" %
1961  B0signalPercentage,
1962  "% B0bars = ",
1963  " %0.4f" %
1964  B0signalbarsPercentage,
1965  "% Diff = ",
1966  " %0.2f" %
1967  Diff,
1968  "%")
1969  print("Total amount of events = ", self.nentries)
1970 
1971 
1972 fitDeltaTModule = fitDeltaT()
1973 validation_path.add_module(fitDeltaTModule)
1974 
1975 b2.process(validation_path)
1976 print(b2.statistics)
B0_GenDeltaTFit.fitDeltaT.DT0
DT0
Definition: B0_GenDeltaTFit.py:335
B0_GenDeltaTFit.fitDeltaT.fitDataTagDeltaTNeg
fitDataTagDeltaTNeg
Definition: B0_GenDeltaTFit.py:354
B0_GenDeltaTFit.fitDeltaT.__init__
def __init__(self)
Definition: B0_GenDeltaTFit.py:324
B0_GenDeltaTFit.fitDeltaT.fitData
fitData
Definition: B0_GenDeltaTFit.py:351
B0_GenDeltaTFit.fitDeltaT.DeltaZtagUpsilon
DeltaZtagUpsilon
Definition: B0_GenDeltaTFit.py:370
B0_GenDeltaTFit.fitDeltaT.fitDataSig
fitDataSig
Definition: B0_GenDeltaTFit.py:355
B0_GenDeltaTFit.fitDeltaT.DT
DT
Definition: B0_GenDeltaTFit.py:332
B0_GenDeltaTFit.fitDeltaT.Tau
Tau
Definition: B0_GenDeltaTFit.py:336
B0_GenDeltaTFit.fitDeltaT.nentries
nentries
Definition: B0_GenDeltaTFit.py:330
B0_GenDeltaTFit.fitDeltaT.B0sInSignalSide
B0sInSignalSide
Definition: B0_GenDeltaTFit.py:362
B0_GenDeltaTFit.fitDeltaT.fitDataDeltaZtagUpsilon
fitDataDeltaZtagUpsilon
Definition: B0_GenDeltaTFit.py:383
Belle2::PyStoreObj
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:69
B0_GenDeltaTFit.fitDeltaT.fitDataTagDeltaTPos
fitDataTagDeltaTPos
Definition: B0_GenDeltaTFit.py:353
B0_GenDeltaTFit.fitDeltaT.A
A
Definition: B0_GenDeltaTFit.py:337
B0_GenDeltaTFit.fitDeltaT.BetaGamma
BetaGamma
Definition: B0_GenDeltaTFit.py:341
B0_GenDeltaTFit.fitDeltaT.B0barsInSignalSide
B0barsInSignalSide
Definition: B0_GenDeltaTFit.py:363
B0_GenDeltaTFit.fitDeltaT.TsigPosNeg
TsigPosNeg
Definition: B0_GenDeltaTFit.py:392
B0_GenDeltaTFit.fitDeltaT.Ttag
Ttag
Definition: B0_GenDeltaTFit.py:334
B0_GenDeltaTFit.fitDeltaT.S
S
Definition: B0_GenDeltaTFit.py:338
B0_GenDeltaTFit.fitDeltaT.fitDataTag
fitDataTag
Definition: B0_GenDeltaTFit.py:352
B0_GenDeltaTFit.fitDeltaT
Definition: B0_GenDeltaTFit.py:252
B0_GenDeltaTFit.fitDeltaT.event
def event(self)
Definition: B0_GenDeltaTFit.py:395
B0_GenDeltaTFit.fitDeltaT.fitDataSigDeltaTPos
fitDataSigDeltaTPos
Definition: B0_GenDeltaTFit.py:356
B0_GenDeltaTFit.fitDeltaT.fractionB0barInSignalSide
fractionB0barInSignalSide
Definition: B0_GenDeltaTFit.py:364
B0_GenDeltaTFit.fitDeltaT.fitDataDeltaZsigUpsilon
fitDataDeltaZsigUpsilon
Definition: B0_GenDeltaTFit.py:377
B0_GenDeltaTFit.fitDeltaT.B0barsInTagSide
B0barsInTagSide
Definition: B0_GenDeltaTFit.py:360
B0_GenDeltaTFit.fitDeltaT.TtagPosNeg
TtagPosNeg
Definition: B0_GenDeltaTFit.py:393
B0_GenDeltaTFit.fitDeltaT.DeltaZsigUpsilon
DeltaZsigUpsilon
Definition: B0_GenDeltaTFit.py:369
B0_GenDeltaTFit.fitDeltaT.Norm
Norm
Definition: B0_GenDeltaTFit.py:340
B0_GenDeltaTFit.fitDeltaT.qsig
qsig
Definition: B0_GenDeltaTFit.py:347
B0_GenDeltaTFit.fitDeltaT.Tsig
Tsig
Definition: B0_GenDeltaTFit.py:333
B0_GenDeltaTFit.fitDeltaT.B0sInTagSide
B0sInTagSide
Definition: B0_GenDeltaTFit.py:359
B0_GenDeltaTFit.fitDeltaT.q
q
Definition: B0_GenDeltaTFit.py:343
B0_GenDeltaTFit.fitDeltaT.terminate
def terminate(self)
Definition: B0_GenDeltaTFit.py:495
B0_GenDeltaTFit.fitDeltaT.fractionB0InSignalSide
fractionB0InSignalSide
Definition: B0_GenDeltaTFit.py:365
B0_GenDeltaTFit.fitDeltaT.fitDataSigDeltaTNeg
fitDataSigDeltaTNeg
Definition: B0_GenDeltaTFit.py:357
B0_GenDeltaTFit.fitDeltaT.DM
DM
Definition: B0_GenDeltaTFit.py:339