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