Belle II Software  light-2403-persian
B0_GenDeltaTFit.py
1 #!/usr/bin/env python
2 
3 
10 
11 
25 
26 import ROOT
27 import basf2 as b2
28 from ROOT import Belle2
29 import modularAnalysis as ma
30 import numpy as np
31 import math
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(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().__init__()
330  self.nentriesnentries = 0
331 
332  self.DTDT = ROOT.RooRealVar("DT", "#Delta#it{t}", 0., -11., 11., "ps")
333  self.TsigTsig = ROOT.RooRealVar("Tsig", "#it{t}_{sig}^{gen}", 0., 0., 11., "ps")
334  self.TtagTtag = ROOT.RooRealVar("Ttag", "#it{t}_{tag}^{gen}", 0., 0., 11., "ps")
335  self.DT0DT0 = ROOT.RooRealVar("DT0", "DT0", 0., -7., 7.)
336  self.TauTau = ROOT.RooRealVar("Tau", "Tau", 1.525, 0., 4.)
337  self.AA = ROOT.RooRealVar("A", "A", 0, -1, 1)
338  self.SS = ROOT.RooRealVar("S", "S", 0, -1, 1)
339  self.DMDM = ROOT.RooRealVar("DM", "DM", 0.507)
340  self.NormNorm = ROOT.RooRealVar("Norm", "Norm", 2, 1, 8)
341  self.BetaGammaBetaGamma = ROOT.RooRealVar("BetaGamma", "BetaGamma", 0.5, 0.1, 1)
342 
343  self.qq = ROOT.RooCategory("q", "q")
344  self.qq.defineType("1")
345  self.qq.defineType("-1")
346 
347  self.qsigqsig = ROOT.RooCategory("qsig", "qsig")
348  self.qsigqsig.defineType("1")
349  self.qsigqsig.defineType("-1")
350 
351  self.fitDatafitData = ROOT.RooDataSet("fitData", "fitData", ROOT.RooArgSet(self.DTDT, self.qq))
352  self.fitDataTagfitDataTag = ROOT.RooDataSet("fitDataTag", "fitDataTag", ROOT.RooArgSet(self.TtagTtag, self.qq))
353  self.fitDataTagDeltaTPosfitDataTagDeltaTPos = ROOT.RooDataSet("fitDataTagDeltaTPos", "fitDataTagDeltaTPos", ROOT.RooArgSet(self.TtagTtag, self.qq))
354  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg = ROOT.RooDataSet("fitDataTagDeltaTNeg", "fitDataTagDeltaTNeg", ROOT.RooArgSet(self.TtagTtag, self.qq))
355  self.fitDataSigfitDataSig = ROOT.RooDataSet("fitDataSig", "fitDataSig", ROOT.RooArgSet(self.TsigTsig, self.qq, self.qsigqsig))
356  self.fitDataSigDeltaTPosfitDataSigDeltaTPos = ROOT.RooDataSet("fitDataSigDeltaTPos", "fitDataSigDeltaTPos", ROOT.RooArgSet(self.TsigTsig, self.qq))
357  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg = ROOT.RooDataSet("fitDataSigDeltaTNeg", "fitDataSigDeltaTNeg", ROOT.RooArgSet(self.TsigTsig, self.qq))
358 
359  self.B0sInTagSideB0sInTagSide = 0
360  self.B0barsInTagSideB0barsInTagSide = 0
361 
362  self.B0sInSignalSideB0sInSignalSide = 0
363  self.B0barsInSignalSideB0barsInSignalSide = 0
364  self.fractionB0barInSignalSidefractionB0barInSignalSide = 0.7057091319609885
365  self.fractionB0InSignalSidefractionB0InSignalSide = 0.2942908680390115
366 
367  # Study od Delta z
368 
369  self.DeltaZsigUpsilonDeltaZsigUpsilon = ROOT.RooRealVar("DeltaZsigUpsilon", "(#it{z}_{sig}^{dec} - #it{z}_{sig}^{prod})^{gen}", 0., 0., 0.8)
370  self.DeltaZtagUpsilonDeltaZtagUpsilon = 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.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon = ROOT.RooDataSet(
378  "fitDataDeltaZsigUpsilon",
379  "fitDataDeltaZsigUpsilon",
380  ROOT.RooArgSet(
381  self.DeltaZsigUpsilonDeltaZsigUpsilon,
382  self.qq))
383  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon = ROOT.RooDataSet(
384  "fitDataDeltaZtagUpsilon",
385  "fitDataDeltaZtagUpsilon",
386  ROOT.RooArgSet(
387  self.DeltaZtagUpsilonDeltaZtagUpsilon,
388  self.qq))
389 
390  # Only for plotting
391 
392  self.TsigPosNegTsigPosNeg = ROOT.RooRealVar("Tsig", "#it{t}_{sig}^{gen}", 0., 0., 7., "ps")
393  self.TtagPosNegTtagPosNeg = 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  rejectEventFlag = False
450 
451  DeltaT = 1000 * (tSig - tTag)
452  # DeltaZ = (zSig - zTag)
453  # DeltaT = 100*DeltaZ /(0.284*3)
454  if len(B0Flags) == 2:
455  if B0Flags[0] != B0Flags[1] and rejectEventFlag is not True:
456  self.nentriesnentries += 1
457  # DeltaT = 1000*(tSig - tTag)
458  # print(DeltaT)
459  self.DTDT.setVal(DeltaT)
460  self.TsigTsig.setVal(float(1000 * tSig))
461  self.TtagTtag.setVal(float(1000 * tTag))
462 
463  if tagPDG > 0:
464  self.qq.setLabel("1")
465  self.B0sInTagSideB0sInTagSide += 1
466  elif tagPDG < 0:
467  self.qq.setLabel("-1")
468  self.B0barsInTagSideB0barsInTagSide += 1
469  if sigPDG > 0:
470  self.qsigqsig.setLabel("1")
471  self.B0sInSignalSideB0sInSignalSide += 1
472  elif sigPDG < 0:
473  self.qsigqsig.setLabel("-1")
474  self.B0barsInSignalSideB0barsInSignalSide += 1
475  self.fitDatafitData.add(ROOT.RooArgSet(self.DTDT, self.qq))
476  self.fitDataSigfitDataSig.add(ROOT.RooArgSet(self.TsigTsig, self.qq, self.qsigqsig))
477  self.fitDataTagfitDataTag.add(ROOT.RooArgSet(self.TtagTtag, self.qq))
478 
479  if DeltaT >= 0:
480  self.fitDataTagDeltaTPosfitDataTagDeltaTPos.add(ROOT.RooArgSet(self.TtagTtag, self.qq))
481  self.fitDataSigDeltaTPosfitDataSigDeltaTPos.add(ROOT.RooArgSet(self.TsigTsig, self.qq))
482 
483  if DeltaT <= 0:
484  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg.add(ROOT.RooArgSet(self.TtagTtag, self.qq))
485  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg.add(ROOT.RooArgSet(self.TsigTsig, self.qq))
486 
487  self.DeltaZtagUpsilonDeltaZtagUpsilon.setVal(10 * float(zTag - zUpsilon))
488  self.DeltaZsigUpsilonDeltaZsigUpsilon.setVal(10 * float(zSig - zUpsilon))
489  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon.add(ROOT.RooArgSet(self.DeltaZtagUpsilonDeltaZtagUpsilon, self.qq))
490  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon.add(ROOT.RooArgSet(self.DeltaZsigUpsilonDeltaZsigUpsilon, self.qq))
491 
492  def terminate(self):
493  """
494  Here the known quantum mechanical pdfs are defined and fitted to the generated MC distributions.
495  Afterwards, the plots are saved.
496  """
497  # S
498 
499  # Delta t Fit
500 
501  self.fitDatafitData.Print()
502 
503  self.DMDM.setConstant(ROOT.kTRUE)
504  self.TauTau.setConstant(ROOT.kTRUE)
505  self.DT0DT0.setConstant(ROOT.kTRUE)
506  self.NormNorm.setConstant(ROOT.kTRUE)
507  # self.A.setConstant(ROOT.kTRUE)
508  # self.S.setConstant(ROOT.kTRUE)
509 
510  # self.S.setVal(6.9184e-01)
511  # self.S.setConstant(ROOT.kTRUE)
512  # self.A.setVal(4.4702e-03)
513  # self.A.setConstant(ROOT.kTRUE)
514 
515  DeltaTp = ROOT.RooGenericPdf(
516  "DeltaTp",
517  "DeltaTp",
518  "(exp(-1*abs(DT-DT0)/Tau)/(2*Norm*Tau))*(1+(A*cos(DM*DT)+S*sin(DM*DT)))",
519  ROOT.RooArgList(
520  self.NormNorm,
521  self.DTDT,
522  self.DT0DT0,
523  self.TauTau,
524  self.AA,
525  self.SS,
526  self.DMDM))
527  DeltaTm = ROOT.RooGenericPdf(
528  "DeltaTm",
529  "DeltaTm",
530  "(exp(-1*abs(DT-DT0)/Tau)/(2*Norm*Tau))*(1-(A*cos(DM*DT)+S*sin(DM*DT)))",
531  ROOT.RooArgList(
532  self.NormNorm,
533  self.DTDT,
534  self.DT0DT0,
535  self.TauTau,
536  self.AA,
537  self.SS,
538  self.DMDM))
539  DeltaT = ROOT.RooSimultaneous("DeltaT", "DeltaT", self.qq)
540  DeltaT.addPdf(DeltaTp, "1")
541  DeltaT.addPdf(DeltaTm, "-1")
542  DeltaT.Print()
543 
544  fitRes = DeltaT.fitTo(
545  self.fitDatafitData, ROOT.RooFit.Minos(
546  ROOT.kFALSE), ROOT.RooFit.Extended(
547  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
548  fitRes.Print("v")
549 
550  DTFrame = self.DTDT.frame()
551  print("here", str(DTFrame.GetXaxis().GetXmax()))
552  self.fitDatafitData.plotOn(DTFrame, ROOT.RooFit.Cut("abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())))
553  self.fitDatafitData.plotOn(DTFrame,
554  ROOT.RooFit.Cut("q==q::1" + " && abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())),
555  ROOT.RooFit.MarkerColor(ROOT.kRed),
556  ROOT.RooFit.Name("data_histB0"))
557  self.fitDatafitData.plotOn(DTFrame, ROOT.RooFit.Cut("q==1" + " && abs(DT) < " + str(DTFrame.GetXaxis().GetXmax())),
558  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1), ROOT.RooFit.Name("data_histAB0"))
559 
560  # DeltaT.plotOn(DTFrame, ROOT.RooFit.Slice(self.q,"1"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q),self.fitData),
561  # ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor( ROOT.kRed),ROOT.RooFit.Components("B02JpsiKsp"));
562  # DeltaT.plotOn(DTFrame, ROOT.RooFit.Slice(self.q,"-1"), ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.q),self.fitData),
563  # ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor( ROOT.kBlue+1), ROOT.RooFit.Components("B02JpsiKsm"));
564  DeltaT.plotOn(
565  DTFrame,
566  ROOT.RooFit.ProjWData(
567  ROOT.RooArgSet(
568  self.qq),
569  self.fitDatafitData),
570  ROOT.RooFit.LineWidth(3),
571  ROOT.RooFit.LineColor(
572  ROOT.kBlack),
573  ROOT.RooFit.LineStyle(3))
574  DeltaT.plotOn(
575  DTFrame,
576  ROOT.RooFit.ProjWData(
577  ROOT.RooArgSet(
578  self.qq),
579  self.fitDatafitData),
580  ROOT.RooFit.LineWidth(3),
581  ROOT.RooFit.LineColor(
582  ROOT.kRed),
583  ROOT.RooFit.Components("DeltaTp"))
584  DeltaT.plotOn(DTFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.qq), self.fitDatafitData), ROOT.RooFit.LineWidth(
585  3), ROOT.RooFit.LineColor(ROOT.kBlue + 1), ROOT.RooFit.LineStyle(7), ROOT.RooFit.Components("DeltaTm"))
586 
587  DTFrame.SetTitle("")
588  DTFrame.GetXaxis().SetTitle("#Delta#it{t}^{gen} [ps]")
589  DTFrame.GetXaxis().SetTitleSize(0.05)
590  DTFrame.GetXaxis().SetLabelSize(0.045)
591  DTFrame.GetYaxis().SetTitleSize(0.05)
592  DTFrame.GetYaxis().SetTitleOffset(1.5)
593  DTFrame.GetYaxis().SetLabelSize(0.045)
594 
595  c1 = ROOT.TCanvas("c1", "c1", 1400, 1100)
596  c1.cd()
597 
598  self.DTDT.Print()
599 
600  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
601  Pad.SetLeftMargin(0.15)
602  Pad.SetBottomMargin(0.15)
603  Pad.Draw()
604  Pad.cd()
605  DTFrame.Draw()
606  DTFrame.Print()
607 
608  curveMC = DTFrame.getCurve("DeltaT_Norm[DT]")
609  curveMCB0 = DTFrame.getCurve("DeltaT_Norm[DT]_Comp[DeltaTp]")
610  curveMCAB0 = DTFrame.getCurve("DeltaT_Norm[DT]_Comp[DeltaTm]")
611  curveMC.SetFillColor(ROOT.kWhite)
612  curveMCB0.SetFillColor(ROOT.kWhite)
613  curveMCAB0.SetFillColor(ROOT.kWhite)
614 
615  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
616  lMC.AddEntry(curveMC, "Both")
617  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
618  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
619 
620  lMC.SetTextSize(0.065)
621  lMC.Draw()
622 
623  textFitRes = "#it{A}_{#it{CP}} = " + \
624  f'{self.A.getValV() if abs(self.A.getValV()) > 0.005 else abs(self.A.getValV()): 1.2f}' + \
625  ", " + "#it{S}_{#it{CP}} = " + \
626  f'{self.S.getValV() if abs(self.S.getValV()) > 0.005 else abs(self.S.getValV()): 1.2f}'
627 
628  lFitRes = ROOT.TLegend(0.215, 0.905, 0.7, 0.98)
629  lFitRes.AddEntry("", textFitRes, " ")
630  lFitRes.SetBorderSize(0)
631  lFitRes.SetTextSize(0.054)
632  lFitRes.Draw()
633 
634  # lFitRes = ROOT.TLegend(0.18, 0.68, 0.40, 0.89)
635  # lFitRes.AddEntry("", "#it{A}_{#it{CP}} = " + '{: 1.2f}'.format(self.A.getValV()), " ")
636  # lFitRes.AddEntry("", "#it{S}_{#it{CP}} = " + '{: 1.2f}'.format(self.S.getValV()), " " )
637  # lFitRes.SetBorderSize(0)
638  # lFitRes.SetTextSize(0.054)
639  # lFitRes.Draw()
640 
641  ROOT.TText(
642  0.1,
643  1,
644  "#it{A}_{#it{CP}} = " +
645  f'{self.A.getValV():1.2f}' +
646  " #pm " +
647  f'{self.A.getError():1.2f}')
648  ROOT.TText(
649  0.1,
650  0.5,
651  "#it{A}_{#it{CP}} = " +
652  f'{self.S.getValV():1.2f}' +
653  " #pm " +
654  f'{self.S.getError():1.2f}')
655 
656  fitRes.Print("v")
657  c1.SaveAs("./" + 'B0_JPsiKs_DeltaTFittedOnGenMC.pdf')
658  c1.Clear()
659 
660  dotsMC = DTFrame.getHist("h_fitData")
661  dotsMCB0 = DTFrame.getHist("data_histB0")
662  dotsMCAB0 = DTFrame.getHist("data_histAB0")
663 
664  plotWithAsymmetry(
665  DTFrame,
666  self.DTDT,
667  dotsMC,
668  dotsMCB0,
669  dotsMCAB0,
670  curveMC,
671  curveMCB0,
672  curveMCAB0,
673  "[ps]",
674  "DeltaTFittedOnGenMC",
675  lMC,
676  textFitRes,
677  "#it{a}_{#it{CP}}(#Delta#it{t})",
678  True)
679 
680  # -----------------
681 
682  # z_tag fit
683 
684  # z_tag fit whole Delta t range
685 
686  # self.S.setVal(0.)
687  # self.A.setVal(0.)
688  self.SS.setConstant(ROOT.kTRUE)
689  self.AA.setConstant(ROOT.kTRUE)
690 
691  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon.Print()
692  DeltaZtagUpsilonp = ROOT.RooGenericPdf(
693  "DeltaZtagUpsilonp",
694  "DeltaZtagUpsilonp",
695  "(exp(-1*abs(DeltaZtagUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 + " +
696  "(A*(cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))) + " +
697  " S*((DM*Tau)*cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma)) - sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))))" +
698  "/(1+(Tau*DM)*(Tau*DM)))",
699  ROOT.RooArgList(
700  self.NormNorm,
701  self.BetaGammaBetaGamma,
702  self.DeltaZtagUpsilonDeltaZtagUpsilon,
703  self.TauTau,
704  self.AA,
705  self.SS,
706  self.DMDM))
707  DeltaZtagUpsilonm = ROOT.RooGenericPdf(
708  "DeltaZtagUpsilonm",
709  "DeltaZtagUpsilonm",
710  "(exp(-1*abs(DeltaZtagUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 - " +
711  "(A*(cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))) + " +
712  " S*((DM*Tau)*cos(DM*DeltaZtagUpsilon/(0.3*BetaGamma)) - sin(DM*DeltaZtagUpsilon/(0.3*BetaGamma))))" +
713  "/(1+(Tau*DM)*(Tau*DM)))",
714  ROOT.RooArgList(
715  self.NormNorm,
716  self.BetaGammaBetaGamma,
717  self.DeltaZtagUpsilonDeltaZtagUpsilon,
718  self.TauTau,
719  self.AA,
720  self.SS,
721  self.DMDM))
722  DeltaZtagUpsilonModel = ROOT.RooSimultaneous("DeltaZtagUpsilonModel", "DeltaZtagUpsilonModel", self.qq)
723  DeltaZtagUpsilonModel.addPdf(DeltaZtagUpsilonp, "1")
724  DeltaZtagUpsilonModel.addPdf(DeltaZtagUpsilonm, "-1")
725  DeltaZtagUpsilonModel.Print()
726 
727  fitResTag = DeltaZtagUpsilonModel.fitTo(
728  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon, ROOT.RooFit.Minos(
729  ROOT.kFALSE), ROOT.RooFit.Extended(
730  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
731  fitResTag.Print("v")
732 
733  DeltaZtagUpsilonFrame = self.DeltaZtagUpsilonDeltaZtagUpsilon.frame()
734  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon.plotOn(DeltaZtagUpsilonFrame)
735  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon.plotOn(
736  DeltaZtagUpsilonFrame,
737  ROOT.RooFit.Cut("q==q::1"),
738  ROOT.RooFit.MarkerColor(
739  ROOT.kRed),
740  ROOT.RooFit.Name("data_histB0tag"))
741  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon.plotOn(
742  DeltaZtagUpsilonFrame,
743  ROOT.RooFit.Cut("q==1"),
744  ROOT.RooFit.MarkerColor(
745  ROOT.kBlue + 1),
746  ROOT.RooFit.Name("data_histAB0tag"))
747 
748  DeltaZtagUpsilonModel.plotOn(
749  DeltaZtagUpsilonFrame,
750  ROOT.RooFit.ProjWData(
751  ROOT.RooArgSet(
752  self.qq),
753  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon),
754  ROOT.RooFit.LineWidth(3),
755  ROOT.RooFit.LineColor(
756  ROOT.kBlack),
757  ROOT.RooFit.LineStyle(3))
758  DeltaZtagUpsilonModel.plotOn(
759  DeltaZtagUpsilonFrame,
760  ROOT.RooFit.ProjWData(
761  ROOT.RooArgSet(
762  self.qq),
763  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon),
764  ROOT.RooFit.LineWidth(3),
765  ROOT.RooFit.LineColor(
766  ROOT.kRed),
767  ROOT.RooFit.Components("DeltaZtagUpsilonp"))
768  DeltaZtagUpsilonModel.plotOn(
769  DeltaZtagUpsilonFrame,
770  ROOT.RooFit.ProjWData(
771  ROOT.RooArgSet(
772  self.qq),
773  self.fitDataDeltaZtagUpsilonfitDataDeltaZtagUpsilon),
774  ROOT.RooFit.LineWidth(3),
775  ROOT.RooFit.LineColor(
776  ROOT.kBlue + 1),
777  ROOT.RooFit.LineStyle(7),
778  ROOT.RooFit.Components("DeltaZtagUpsilonm"))
779  DeltaZtagUpsilonFrame.SetTitle("")
780  DeltaZtagUpsilonFrame.GetXaxis().SetTitle("(#it{z}_{tag}^{dec} - #it{z}_{tag}^{prod}){}^{gen} [mm]")
781  DeltaZtagUpsilonFrame.GetXaxis().SetTitleSize(0.05)
782  DeltaZtagUpsilonFrame.GetXaxis().SetLabelSize(0.045)
783  DeltaZtagUpsilonFrame.GetYaxis().SetTitleSize(0.05)
784  DeltaZtagUpsilonFrame.GetYaxis().SetTitleOffset(1.5)
785  DeltaZtagUpsilonFrame.GetYaxis().SetLabelSize(0.045)
786 
787  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
788  c2.cd()
789 
790  self.DeltaZtagUpsilonDeltaZtagUpsilon.Print()
791 
792  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
793  Pad.SetLeftMargin(0.15)
794  Pad.SetBottomMargin(0.15)
795  Pad.Draw()
796  Pad.cd()
797  DeltaZtagUpsilonFrame.Draw()
798 
799  curveMC = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]")
800  curveMCB0 = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]_Comp[DeltaZtagUpsilonp]")
801  curveMCAB0 = DeltaZtagUpsilonFrame.getCurve("DeltaZtagUpsilonModel_Norm[DeltaZtagUpsilon]_Comp[DeltaZtagUpsilonm]")
802  curveMC.SetFillColor(ROOT.kWhite)
803  curveMCB0.SetFillColor(ROOT.kWhite)
804  curveMCAB0.SetFillColor(ROOT.kWhite)
805 
806  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
807  lMC.AddEntry(curveMC, "Both")
808  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
809  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
810 
811  lMC.SetTextSize(0.065)
812  lMC.Draw()
813 
814  fitResTag.Print("v")
815  c2.SaveAs("./" + 'B0_JPsiKs_DeltaZtagUpsilonFittedOnGenMC.pdf')
816  c2.Clear()
817 
818  dotsMC = DeltaZtagUpsilonFrame.getHist("h_fitDataDeltaZtagUpsilon")
819  dotsMCB0 = DeltaZtagUpsilonFrame.getHist("data_histB0tag")
820  dotsMCAB0 = DeltaZtagUpsilonFrame.getHist("data_histAB0tag")
821 
822  plotWithAsymmetry(
823  DeltaZtagUpsilonFrame,
824  self.DeltaZtagUpsilonDeltaZtagUpsilon,
825  dotsMC,
826  dotsMCB0,
827  dotsMCAB0,
828  curveMC,
829  curveMCB0,
830  curveMCAB0,
831  "[mm]",
832  "DeltaZtagUpsilonFittedOnGenMC",
833  lMC,
834  "",
835  labelAsymmetry,
836  True)
837 
838  self.SS.setConstant(ROOT.kFALSE)
839  self.AA.setConstant(ROOT.kFALSE)
840 
841  # -----------------
842 
843  # t_tag fit
844 
845  # t_tag fit whole Delta t range
846 
847  # self.S.setConstant(ROOT.kTRUE)
848  self.SS.setVal(0.)
849  self.AA.setVal(0.)
850  # self.A.setConstant(ROOT.kTRUE)
851 
852  self.fitDataTagfitDataTag.Print()
853  Ttagp = ROOT.RooGenericPdf(
854  "Ttagp",
855  "Ttagp",
856  "(exp(-1*abs(Ttag)/Tau)/(Norm*Tau))*(1 + " +
857  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM)))",
858  ROOT.RooArgList(
859  self.NormNorm,
860  self.TtagTtag,
861  self.TauTau,
862  self.AA,
863  self.SS,
864  self.DMDM))
865  Ttagm = ROOT.RooGenericPdf(
866  "Ttagm",
867  "Ttagm",
868  "(exp(-1*abs(Ttag)/Tau)/(Norm*Tau))*(1 - " +
869  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM)))",
870  ROOT.RooArgList(
871  self.NormNorm,
872  self.TtagTtag,
873  self.TauTau,
874  self.AA,
875  self.SS,
876  self.DMDM))
877  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.qq)
878  TtagModel.addPdf(Ttagp, "1")
879  TtagModel.addPdf(Ttagm, "-1")
880  TtagModel.Print()
881 
882  fitResTag = TtagModel.fitTo(
883  self.fitDataTagfitDataTag, ROOT.RooFit.Minos(
884  ROOT.kFALSE), ROOT.RooFit.Extended(
885  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
886  fitResTag.Print("v")
887 
888  TtagFrame = self.TtagTtag.frame()
889  self.fitDataTagfitDataTag.plotOn(TtagFrame)
890  self.fitDataTagfitDataTag.plotOn(
891  TtagFrame,
892  ROOT.RooFit.Cut("q==q::1"),
893  ROOT.RooFit.MarkerColor(
894  ROOT.kRed),
895  ROOT.RooFit.Name("data_histB0tag"))
896  self.fitDataTagfitDataTag.plotOn(
897  TtagFrame,
898  ROOT.RooFit.Cut("q==1"),
899  ROOT.RooFit.MarkerColor(
900  ROOT.kBlue + 1),
901  ROOT.RooFit.Name("data_histAB0tag"))
902 
903  TtagModel.plotOn(
904  TtagFrame,
905  ROOT.RooFit.ProjWData(
906  ROOT.RooArgSet(
907  self.qq),
908  self.fitDataTagfitDataTag),
909  ROOT.RooFit.LineWidth(3),
910  ROOT.RooFit.LineColor(
911  ROOT.kBlack),
912  ROOT.RooFit.LineStyle(3))
913  TtagModel.plotOn(
914  TtagFrame,
915  ROOT.RooFit.ProjWData(
916  ROOT.RooArgSet(
917  self.qq),
918  self.fitDataTagfitDataTag),
919  ROOT.RooFit.LineWidth(3),
920  ROOT.RooFit.LineColor(
921  ROOT.kRed),
922  ROOT.RooFit.Components("Ttagp"))
923  TtagModel.plotOn(
924  TtagFrame,
925  ROOT.RooFit.ProjWData(
926  ROOT.RooArgSet(
927  self.qq),
928  self.fitDataTagfitDataTag),
929  ROOT.RooFit.LineWidth(3),
930  ROOT.RooFit.LineColor(
931  ROOT.kBlue + 1),
932  ROOT.RooFit.LineStyle(7),
933  ROOT.RooFit.Components("Ttagm"))
934  TtagFrame.SetTitle("")
935  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
936  TtagFrame.GetXaxis().SetTitleSize(0.05)
937  TtagFrame.GetXaxis().SetLabelSize(0.045)
938  TtagFrame.GetYaxis().SetTitleSize(0.05)
939  TtagFrame.GetYaxis().SetTitleOffset(1.5)
940  TtagFrame.GetYaxis().SetLabelSize(0.045)
941 
942  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
943  c2.cd()
944 
945  self.TtagTtag.Print()
946 
947  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
948  Pad.SetLeftMargin(0.15)
949  Pad.SetBottomMargin(0.15)
950  Pad.Draw()
951  Pad.cd()
952  TtagFrame.Draw()
953 
954  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
955  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
956  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
957  curveMC.SetFillColor(ROOT.kWhite)
958  curveMCB0.SetFillColor(ROOT.kWhite)
959  curveMCAB0.SetFillColor(ROOT.kWhite)
960 
961  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
962  lMC.AddEntry(curveMC, "Both")
963  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
964  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
965 
966  lMC.SetTextSize(0.065)
967  lMC.Draw()
968 
969  fitResTag.Print("v")
970  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMC.pdf')
971  c2.Clear()
972 
973  dotsMC = TtagFrame.getHist("h_fitDataTag")
974  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
975  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
976 
977  plotWithAsymmetry(TtagFrame, self.TtagTtag, dotsMC, dotsMCB0, dotsMCAB0, curveMC, curveMCB0,
978  curveMCAB0, "[ps]", "TtagFittedOnGenMC", lMC, "", labelAsymmetry, True)
979 
980  # ---
981 
982  # t_tag fit positive Delta t range
983 
984  # self.S.setConstant(ROOT.kTRUE)
985  # self.S.setVal(0.)
986  # self.A.setVal(0.)
987  # self.A.setConstant(ROOT.kTRUE)
988 
989  self.fitDataTagDeltaTPosfitDataTagDeltaTPos.Print()
990  Ttagp = ROOT.RooGenericPdf(
991  "Ttagp",
992  "Ttagp",
993  "(exp(-2*abs(Ttag)/Tau)/(Tau))*(1 + (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
994  ROOT.RooArgList(
995  self.NormNorm,
996  self.TtagTtag,
997  self.TauTau,
998  self.AA,
999  self.SS,
1000  self.DMDM))
1001  Ttagm = ROOT.RooGenericPdf(
1002  "Ttagm",
1003  "Ttagm",
1004  "(exp(-2*abs(Ttag)/Tau)/(Tau))*(1 - (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1005  ROOT.RooArgList(
1006  self.NormNorm,
1007  self.TtagTtag,
1008  self.TauTau,
1009  self.AA,
1010  self.SS,
1011  self.DMDM))
1012  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.qq)
1013  TtagModel.addPdf(Ttagp, "1")
1014  TtagModel.addPdf(Ttagm, "-1")
1015  TtagModel.Print()
1016 
1017  fitResTagPos = TtagModel.fitTo(
1018  self.fitDataTagDeltaTPosfitDataTagDeltaTPos, ROOT.RooFit.Minos(
1019  ROOT.kFALSE), ROOT.RooFit.Extended(
1020  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1021  fitResTagPos.Print("v")
1022 
1023  TtagFrame = self.TtagPosNegTtagPosNeg.frame()
1024  self.fitDataTagDeltaTPosfitDataTagDeltaTPos.plotOn(TtagFrame, ROOT.RooFit.Cut("Ttag < " + str(TtagFrame.GetXaxis().GetXmax())))
1025  self.fitDataTagDeltaTPosfitDataTagDeltaTPos.plotOn(TtagFrame,
1026  ROOT.RooFit.Cut("q==q::1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1027  ROOT.RooFit.MarkerColor(ROOT.kRed),
1028  ROOT.RooFit.Name("data_histB0tag"))
1029  self.fitDataTagDeltaTPosfitDataTagDeltaTPos.plotOn(TtagFrame,
1030  ROOT.RooFit.Cut("q==1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1031  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1032  ROOT.RooFit.Name("data_histAB0tag"))
1033 
1034  TtagModel.plotOn(TtagFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.qq), self.fitDataTagDeltaTPosfitDataTagDeltaTPos),
1035  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1036  TtagModel.plotOn(
1037  TtagFrame,
1038  ROOT.RooFit.ProjWData(
1039  ROOT.RooArgSet(
1040  self.qq),
1041  self.fitDataTagDeltaTPosfitDataTagDeltaTPos),
1042  ROOT.RooFit.LineWidth(3),
1043  ROOT.RooFit.LineColor(
1044  ROOT.kRed),
1045  ROOT.RooFit.Components("Ttagp"))
1046  TtagModel.plotOn(
1047  TtagFrame,
1048  ROOT.RooFit.ProjWData(
1049  ROOT.RooArgSet(
1050  self.qq),
1051  self.fitDataTagDeltaTPosfitDataTagDeltaTPos),
1052  ROOT.RooFit.LineWidth(3),
1053  ROOT.RooFit.LineColor(
1054  ROOT.kBlue + 1),
1055  ROOT.RooFit.LineStyle(7),
1056  ROOT.RooFit.Components("Ttagm"))
1057  TtagFrame.SetTitle("")
1058  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
1059  TtagFrame.GetXaxis().SetTitleSize(0.05)
1060  TtagFrame.GetXaxis().SetLabelSize(0.045)
1061  TtagFrame.GetYaxis().SetTitleSize(0.05)
1062  TtagFrame.GetYaxis().SetTitleOffset(1.5)
1063  TtagFrame.GetYaxis().SetLabelSize(0.045)
1064 
1065  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1066  c2.cd()
1067 
1068  self.TtagTtag.Print()
1069 
1070  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1071  Pad.SetLeftMargin(0.15)
1072  Pad.SetBottomMargin(0.15)
1073  Pad.Draw()
1074  Pad.cd()
1075  TtagFrame.Draw()
1076 
1077  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
1078  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
1079  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
1080  curveMC.SetFillColor(ROOT.kWhite)
1081  curveMCB0.SetFillColor(ROOT.kWhite)
1082  curveMCAB0.SetFillColor(ROOT.kWhite)
1083 
1084  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1085  lMC.AddEntry(curveMC, "Both")
1086  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1087  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1088 
1089  lMC.SetTextSize(0.065)
1090  lMC.Draw()
1091 
1092  fitResTagPos.Print("v")
1093  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMCPositiveDeltaT.pdf')
1094  c2.Clear()
1095 
1096  dotsMC = TtagFrame.getHist("h_fitDataTag")
1097  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
1098  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
1099 
1100  plotWithAsymmetry(
1101  TtagFrame,
1102  self.TtagPosNegTtagPosNeg,
1103  dotsMC,
1104  dotsMCB0,
1105  dotsMCAB0,
1106  curveMC,
1107  curveMCB0,
1108  curveMCAB0,
1109  "[ps]",
1110  "TtagFittedOnGenMCPositiveDeltaT",
1111  lMC,
1112  "",
1113  labelAsymmetry,
1114  True)
1115 
1116  # t_tag fit negative Delta t range
1117 
1118  # self.S.setConstant(ROOT.kTRUE)
1119  self.SS.setVal(0.)
1120  self.AA.setVal(0.)
1121  # self.A.setConstant(ROOT.kTRUE)
1122  self.AA.setConstant(ROOT.kFALSE)
1123 
1124  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg.Print()
1125  Ttagp = ROOT.RooGenericPdf(
1126  "Ttagp",
1127  "Ttagp",
1128  "(exp(-1*abs(Ttag)/Tau)*(1 + " +
1129  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM))) - " +
1130  "(exp(-2*abs(Ttag)/Tau))*(1 + (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1131  ROOT.RooArgList(
1132  self.TtagTtag,
1133  self.TauTau,
1134  self.AA,
1135  self.SS,
1136  self.DMDM))
1137  Ttagm = ROOT.RooGenericPdf(
1138  "Ttagm",
1139  "Ttagm",
1140  "(exp(-1*abs(Ttag)/Tau)*(1 - " +
1141  "(A*(cos(DM*Ttag)+(DM*Tau)*sin(DM*Ttag)) + S*((DM*Tau)*cos(DM*Ttag) - sin(DM*Ttag)))/(1+(Tau*DM)*(Tau*DM))) - " +
1142  "(exp(-2*abs(Ttag)/Tau))*(1 - (A + S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1143  ROOT.RooArgList(
1144  self.TtagTtag,
1145  self.TauTau,
1146  self.AA,
1147  self.SS,
1148  self.DMDM))
1149  TtagModel = ROOT.RooSimultaneous("TtagModel", "TtagModel", self.qq)
1150  TtagModel.addPdf(Ttagp, "1")
1151  TtagModel.addPdf(Ttagm, "-1")
1152  TtagModel.Print()
1153 
1154  fitResTagNeg = TtagModel.fitTo(
1155  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg, ROOT.RooFit.Minos(
1156  ROOT.kFALSE), ROOT.RooFit.Extended(
1157  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1158  fitResTagNeg.Print("v")
1159 
1160  TtagFrame = self.TtagPosNegTtagPosNeg.frame()
1161  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg.plotOn(TtagFrame, ROOT.RooFit.Cut("Ttag < " + str(TtagFrame.GetXaxis().GetXmax())))
1162  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg.plotOn(TtagFrame,
1163  ROOT.RooFit.Cut("q==q::1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1164  ROOT.RooFit.MarkerColor(ROOT.kRed),
1165  ROOT.RooFit.Name("data_histB0tag"))
1166  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg.plotOn(TtagFrame,
1167  ROOT.RooFit.Cut("q==1 && Ttag < " + str(TtagFrame.GetXaxis().GetXmax())),
1168  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1169  ROOT.RooFit.Name("data_histAB0tag"))
1170 
1171  TtagModel.plotOn(TtagFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.qq), self.fitDataTagDeltaTNegfitDataTagDeltaTNeg),
1172  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1173  TtagModel.plotOn(
1174  TtagFrame,
1175  ROOT.RooFit.ProjWData(
1176  ROOT.RooArgSet(
1177  self.qq),
1178  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg),
1179  ROOT.RooFit.LineWidth(3),
1180  ROOT.RooFit.LineColor(
1181  ROOT.kRed),
1182  ROOT.RooFit.Components("Ttagp"))
1183  TtagModel.plotOn(
1184  TtagFrame,
1185  ROOT.RooFit.ProjWData(
1186  ROOT.RooArgSet(
1187  self.qq),
1188  self.fitDataTagDeltaTNegfitDataTagDeltaTNeg),
1189  ROOT.RooFit.LineWidth(3),
1190  ROOT.RooFit.LineColor(
1191  ROOT.kBlue + 1),
1192  ROOT.RooFit.LineStyle(7),
1193  ROOT.RooFit.Components("Ttagm"))
1194  TtagFrame.SetTitle("")
1195  TtagFrame.GetXaxis().SetTitle("#it{t}_{tag}^{gen} [ps]")
1196  TtagFrame.GetXaxis().SetTitleSize(0.05)
1197  TtagFrame.GetXaxis().SetLabelSize(0.045)
1198  TtagFrame.GetYaxis().SetTitleSize(0.05)
1199  TtagFrame.GetYaxis().SetTitleOffset(1.5)
1200  TtagFrame.GetYaxis().SetLabelSize(0.045)
1201 
1202  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1203  c2.cd()
1204 
1205  self.TtagTtag.Print()
1206 
1207  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1208  Pad.SetLeftMargin(0.15)
1209  Pad.SetBottomMargin(0.15)
1210  Pad.Draw()
1211  Pad.cd()
1212  TtagFrame.Draw()
1213  TtagFrame.Print()
1214 
1215  # lMC = ROOT.TLegend(0.7, 0.63, 0.9, 0.89)
1216 
1217  curveMC = TtagFrame.getCurve("TtagModel_Norm[Ttag]")
1218  curveMCB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagp]")
1219  curveMCAB0 = TtagFrame.getCurve("TtagModel_Norm[Ttag]_Comp[Ttagm]")
1220  curveMC.SetFillColor(ROOT.kWhite)
1221  curveMCB0.SetFillColor(ROOT.kWhite)
1222  curveMCAB0.SetFillColor(ROOT.kWhite)
1223 
1224  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1225  lMC.AddEntry(curveMC, "Both")
1226  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1227  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1228 
1229  lMC.SetTextSize(0.065)
1230  lMC.Draw()
1231 
1232  fitResTagNeg.Print("v")
1233  c2.SaveAs("./" + 'B0_JPsiKs_TtagFittedOnGenMCNegativeDeltaT.pdf')
1234  c2.Clear()
1235 
1236  dotsMC = TtagFrame.getHist("h_fitDataTag")
1237  dotsMCB0 = TtagFrame.getHist("data_histB0tag")
1238  dotsMCAB0 = TtagFrame.getHist("data_histAB0tag")
1239 
1240  plotWithAsymmetry(
1241  TtagFrame,
1242  self.TtagPosNegTtagPosNeg,
1243  dotsMC,
1244  dotsMCB0,
1245  dotsMCAB0,
1246  curveMC,
1247  curveMCB0,
1248  curveMCAB0,
1249  "[ps]",
1250  "TtagFittedOnGenMCNegativeDeltaT",
1251  lMC,
1252  "",
1253  labelAsymmetry,
1254  True)
1255 
1256  # ----------------------
1257 
1258  # z_sig fit
1259 
1260  # z_sig fit whole Delta t range
1261 
1262  # self.S.setVal(0.)
1263  # self.A.setVal(0.)
1264  self.SS.setConstant(ROOT.kTRUE)
1265  self.AA.setConstant(ROOT.kTRUE)
1266 
1267  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon.Print()
1268  DeltaZsigUpsilonp = ROOT.RooGenericPdf(
1269  "DeltaZsigUpsilonp",
1270  "DeltaZsigUpsilonp",
1271  "(exp(-1*abs(DeltaZsigUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 + " +
1272  "(A*(cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))) + " +
1273  "S*(-(DM*Tau)*cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma)) + sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))))" +
1274  "/(1+(Tau*DM)*(Tau*DM)))",
1275  ROOT.RooArgList(
1276  self.NormNorm,
1277  self.BetaGammaBetaGamma,
1278  self.DeltaZsigUpsilonDeltaZsigUpsilon,
1279  self.TauTau,
1280  self.AA,
1281  self.SS,
1282  self.DMDM))
1283  DeltaZsigUpsilonm = ROOT.RooGenericPdf(
1284  "DeltaZsigUpsilonm",
1285  "DeltaZsigUpsilonm",
1286  "(exp(-1*abs(DeltaZsigUpsilon/(0.3*BetaGamma))/Tau)/(Norm*Tau))*(1 - "
1287  "(A*(cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma))+(DM*Tau)*sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))) + "
1288  "S*(-(DM*Tau)*cos(DM*DeltaZsigUpsilon/(0.3*BetaGamma)) + sin(DM*DeltaZsigUpsilon/(0.3*BetaGamma))))"
1289  "/(1+(Tau*DM)*(Tau*DM)))",
1290  ROOT.RooArgList(
1291  self.NormNorm,
1292  self.BetaGammaBetaGamma,
1293  self.DeltaZsigUpsilonDeltaZsigUpsilon,
1294  self.TauTau,
1295  self.AA,
1296  self.SS,
1297  self.DMDM))
1298  DeltaZsigUpsilonModel = ROOT.RooSimultaneous("DeltaZsigUpsilonModel", "DeltaZsigUpsilonModel", self.qq)
1299  DeltaZsigUpsilonModel.addPdf(DeltaZsigUpsilonp, "1")
1300  DeltaZsigUpsilonModel.addPdf(DeltaZsigUpsilonm, "-1")
1301  DeltaZsigUpsilonModel.Print()
1302 
1303  fitResTag = DeltaZsigUpsilonModel.fitTo(
1304  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon, ROOT.RooFit.Minos(
1305  ROOT.kFALSE), ROOT.RooFit.Extended(
1306  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1307  fitResTag.Print("v")
1308 
1309  DeltaZsigUpsilonFrame = self.DeltaZsigUpsilonDeltaZsigUpsilon.frame()
1310  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon.plotOn(DeltaZsigUpsilonFrame)
1311  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon.plotOn(
1312  DeltaZsigUpsilonFrame,
1313  ROOT.RooFit.Cut("q==q::1"),
1314  ROOT.RooFit.MarkerColor(
1315  ROOT.kRed),
1316  ROOT.RooFit.Name("data_histB0sig"))
1317  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon.plotOn(
1318  DeltaZsigUpsilonFrame,
1319  ROOT.RooFit.Cut("q==1"),
1320  ROOT.RooFit.MarkerColor(
1321  ROOT.kBlue + 1),
1322  ROOT.RooFit.Name("data_histAB0sig"))
1323 
1324  DeltaZsigUpsilonModel.plotOn(
1325  DeltaZsigUpsilonFrame,
1326  ROOT.RooFit.ProjWData(
1327  ROOT.RooArgSet(
1328  self.qq),
1329  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon),
1330  ROOT.RooFit.LineWidth(3),
1331  ROOT.RooFit.LineColor(
1332  ROOT.kBlack),
1333  ROOT.RooFit.LineStyle(3))
1334  DeltaZsigUpsilonModel.plotOn(
1335  DeltaZsigUpsilonFrame,
1336  ROOT.RooFit.ProjWData(
1337  ROOT.RooArgSet(
1338  self.qq),
1339  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon),
1340  ROOT.RooFit.LineWidth(3),
1341  ROOT.RooFit.LineColor(
1342  ROOT.kRed),
1343  ROOT.RooFit.Components("DeltaZsigUpsilonp"))
1344  DeltaZsigUpsilonModel.plotOn(
1345  DeltaZsigUpsilonFrame,
1346  ROOT.RooFit.ProjWData(
1347  ROOT.RooArgSet(
1348  self.qq),
1349  self.fitDataDeltaZsigUpsilonfitDataDeltaZsigUpsilon),
1350  ROOT.RooFit.LineWidth(3),
1351  ROOT.RooFit.LineColor(
1352  ROOT.kBlue + 1),
1353  ROOT.RooFit.LineStyle(7),
1354  ROOT.RooFit.Components("DeltaZsigUpsilonm"))
1355  DeltaZsigUpsilonFrame.SetTitle("")
1356  DeltaZsigUpsilonFrame.GetXaxis().SetTitle("(#it{z}_{sig}^{dec} - #it{z}_{sig}^{prod}){}^{gen} [mm]")
1357  DeltaZsigUpsilonFrame.GetXaxis().SetTitleSize(0.05)
1358  DeltaZsigUpsilonFrame.GetXaxis().SetLabelSize(0.045)
1359  DeltaZsigUpsilonFrame.GetYaxis().SetTitleSize(0.05)
1360  DeltaZsigUpsilonFrame.GetYaxis().SetTitleOffset(1.5)
1361  DeltaZsigUpsilonFrame.GetYaxis().SetLabelSize(0.045)
1362 
1363  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1364  c2.cd()
1365 
1366  self.DeltaZsigUpsilonDeltaZsigUpsilon.Print()
1367 
1368  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1369  Pad.SetLeftMargin(0.15)
1370  Pad.SetBottomMargin(0.15)
1371  Pad.Draw()
1372  Pad.cd()
1373  DeltaZsigUpsilonFrame.Draw()
1374 
1375  curveMC = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]")
1376  curveMCB0 = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]_Comp[DeltaZsigUpsilonp]")
1377  curveMCAB0 = DeltaZsigUpsilonFrame.getCurve("DeltaZsigUpsilonModel_Norm[DeltaZsigUpsilon]_Comp[DeltaZsigUpsilonm]")
1378  curveMC.SetFillColor(ROOT.kWhite)
1379  curveMCB0.SetFillColor(ROOT.kWhite)
1380  curveMCAB0.SetFillColor(ROOT.kWhite)
1381 
1382  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1383  lMC.AddEntry(curveMC, "Both")
1384  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1385  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1386 
1387  lMC.SetTextSize(0.065)
1388  lMC.Draw()
1389 
1390  fitResTag.Print("v")
1391  c2.SaveAs("./" + 'B0_JPsiKs_DeltaZsigUpsilonFittedOnGenMC.pdf')
1392  c2.Clear()
1393 
1394  dotsMC = DeltaZsigUpsilonFrame.getHist("h_fitDataDeltaZsigUpsilon")
1395  dotsMCB0 = DeltaZsigUpsilonFrame.getHist("data_histB0sig")
1396  dotsMCAB0 = DeltaZsigUpsilonFrame.getHist("data_histAB0sig")
1397 
1398  plotWithAsymmetry(
1399  DeltaZsigUpsilonFrame,
1400  self.DeltaZsigUpsilonDeltaZsigUpsilon,
1401  dotsMC,
1402  dotsMCB0,
1403  dotsMCAB0,
1404  curveMC,
1405  curveMCB0,
1406  curveMCAB0,
1407  "[mm]",
1408  "DeltaZsigUpsilonFittedOnGenMC",
1409  lMC,
1410  "",
1411  labelAsymmetry,
1412  True)
1413 
1414  self.SS.setConstant(ROOT.kFALSE)
1415  self.AA.setConstant(ROOT.kFALSE)
1416 
1417  # --- Signal time fit --------
1418 
1419  self.fitDataSigfitDataSig.Print()
1420  Tsigp = ROOT.RooGenericPdf(
1421  "Tsigp",
1422  "Tsigp",
1423  "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))*(1 + " +
1424  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM)))",
1425  ROOT.RooArgList(
1426  self.NormNorm,
1427  self.TsigTsig,
1428  self.TauTau,
1429  self.AA,
1430  self.SS,
1431  self.DMDM))
1432  Tsigm = ROOT.RooGenericPdf(
1433  "Tsigm",
1434  "Tsigm",
1435  "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))*(1 - " +
1436  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM)))",
1437  ROOT.RooArgList(
1438  self.NormNorm,
1439  self.TsigTsig,
1440  self.TauTau,
1441  self.AA,
1442  self.SS,
1443  self.DMDM))
1444  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.qq)
1445  TsigModel.addPdf(Tsigp, "1")
1446  TsigModel.addPdf(Tsigm, "-1")
1447  TsigModel.Print()
1448 
1449  fitResSig = TsigModel.fitTo(
1450  self.fitDataSigfitDataSig, ROOT.RooFit.Minos(
1451  ROOT.kFALSE), ROOT.RooFit.Extended(
1452  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1453  fitResSig.Print("v")
1454 
1455  TsigFrame = self.TsigTsig.frame()
1456  self.fitDataSigfitDataSig.plotOn(TsigFrame)
1457  self.fitDataSigfitDataSig.plotOn(
1458  TsigFrame,
1459  ROOT.RooFit.Cut("q==q::1"),
1460  ROOT.RooFit.MarkerColor(
1461  ROOT.kRed),
1462  ROOT.RooFit.Name("data_histB0sig"))
1463  self.fitDataSigfitDataSig.plotOn(
1464  TsigFrame,
1465  ROOT.RooFit.Cut("q==1"),
1466  ROOT.RooFit.MarkerColor(
1467  ROOT.kBlue + 1),
1468  ROOT.RooFit.Name("data_histAB0sig"))
1469 
1470  TsigModel.plotOn(
1471  TsigFrame,
1472  ROOT.RooFit.ProjWData(
1473  ROOT.RooArgSet(
1474  self.qsigqsig),
1475  self.fitDataSigfitDataSig),
1476  ROOT.RooFit.LineWidth(3),
1477  ROOT.RooFit.LineColor(
1478  ROOT.kBlack),
1479  ROOT.RooFit.LineStyle(3))
1480  TsigModel.plotOn(
1481  TsigFrame,
1482  ROOT.RooFit.ProjWData(
1483  ROOT.RooArgSet(
1484  self.qsigqsig),
1485  self.fitDataSigfitDataSig),
1486  ROOT.RooFit.LineWidth(3),
1487  ROOT.RooFit.LineColor(
1488  ROOT.kRed),
1489  ROOT.RooFit.Components("Tsigp"))
1490  TsigModel.plotOn(
1491  TsigFrame,
1492  ROOT.RooFit.ProjWData(
1493  ROOT.RooArgSet(
1494  self.qsigqsig),
1495  self.fitDataSigfitDataSig),
1496  ROOT.RooFit.LineWidth(3),
1497  ROOT.RooFit.LineColor(
1498  ROOT.kBlue + 1),
1499  ROOT.RooFit.LineStyle(7),
1500  ROOT.RooFit.Components("Tsigm"))
1501  TsigFrame.SetTitle("")
1502  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1503  TsigFrame.GetXaxis().SetTitleSize(0.05)
1504  TsigFrame.GetXaxis().SetLabelSize(0.045)
1505  TsigFrame.GetYaxis().SetTitleSize(0.05)
1506  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1507  TsigFrame.GetYaxis().SetLabelSize(0.045)
1508 
1509  c3 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1510  c3.cd()
1511 
1512  self.TsigTsig.Print()
1513 
1514  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1515  Pad.SetLeftMargin(0.15)
1516  Pad.SetBottomMargin(0.15)
1517  Pad.Draw()
1518  Pad.cd()
1519  TsigFrame.Draw()
1520 
1521  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1522  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1523  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1524  curveMC.SetFillColor(ROOT.kWhite)
1525  curveMCB0.SetFillColor(ROOT.kWhite)
1526  curveMCAB0.SetFillColor(ROOT.kWhite)
1527 
1528  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1529  lMC.AddEntry(curveMC, "Both")
1530  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1531  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1532 
1533  lMC.SetTextSize(0.065)
1534  lMC.Draw()
1535 
1536  fitResSig.Print("v")
1537  c3.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMC.pdf')
1538  c3.Clear()
1539 
1540  dotsMC = TsigFrame.getHist("h_fitDataSig")
1541  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1542  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1543 
1544  plotWithAsymmetry(TsigFrame, self.TsigTsig, dotsMC, dotsMCB0, dotsMCAB0, curveMC, curveMCB0,
1545  curveMCAB0, "[ps]", "TsigFittedOnGenMC", lMC, "", labelAsymmetry, True)
1546 
1547  # --- ------
1548 
1549  # t_sig fit positive Delta t range
1550 
1551  # self.S.setConstant(ROOT.kTRUE)
1552  # self.S.setVal(0.)
1553  # self.A.setVal(0.)
1554  # self.A.setConstant(ROOT.kTRUE)
1555 
1556  self.fitDataSigDeltaTPosfitDataSigDeltaTPos.Print()
1557 
1558  Tsigp = ROOT.RooGenericPdf(
1559  "Tsigp",
1560  "Tsigp",
1561  "(exp(-1*abs(Tsig)/Tau)*(1 + " +
1562  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM))) - " +
1563  "(exp(-2*abs(Tsig)/Tau))*(1 + (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1564  ROOT.RooArgList(
1565  self.TsigTsig,
1566  self.TauTau,
1567  self.AA,
1568  self.SS,
1569  self.DMDM))
1570  Tsigm = ROOT.RooGenericPdf(
1571  "Tsigm",
1572  "Tsigm",
1573  "(exp(-1*abs(Tsig)/Tau)*(1 - " +
1574  "(A*(cos(DM*Tsig)+(DM*Tau)*sin(DM*Tsig)) + S*(-(DM*Tau)*cos(DM*Tsig) + sin(DM*Tsig)))/(1+(Tau*DM)*(Tau*DM))) - " +
1575  "(exp(-2*abs(Tsig)/Tau))*(1 - (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM))))/(Tau)",
1576  ROOT.RooArgList(
1577  self.TsigTsig,
1578  self.TauTau,
1579  self.AA,
1580  self.SS,
1581  self.DMDM))
1582 
1583  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.qq)
1584  TsigModel.addPdf(Tsigp, "1")
1585  TsigModel.addPdf(Tsigm, "-1")
1586  TsigModel.Print()
1587 
1588  fitResSigPos = TsigModel.fitTo(
1589  self.fitDataSigDeltaTPosfitDataSigDeltaTPos, ROOT.RooFit.Minos(
1590  ROOT.kFALSE), ROOT.RooFit.Extended(
1591  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1592  fitResSigPos.Print("v")
1593 
1594  TsigFrame = self.TsigPosNegTsigPosNeg.frame()
1595  self.fitDataSigDeltaTPosfitDataSigDeltaTPos.plotOn(TsigFrame, ROOT.RooFit.Cut("Tsig < " + str(TsigFrame.GetXaxis().GetXmax())))
1596  self.fitDataSigDeltaTPosfitDataSigDeltaTPos.plotOn(TsigFrame,
1597  ROOT.RooFit.Cut("q==q::1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1598  ROOT.RooFit.MarkerColor(ROOT.kRed),
1599  ROOT.RooFit.Name("data_histB0sig"))
1600  self.fitDataSigDeltaTPosfitDataSigDeltaTPos.plotOn(TsigFrame,
1601  ROOT.RooFit.Cut("q==1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1602  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1603  ROOT.RooFit.Name("data_histAB0sig"))
1604 
1605  TsigModel.plotOn(TsigFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.qq), self.fitDataSigDeltaTPosfitDataSigDeltaTPos),
1606  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1607  TsigModel.plotOn(
1608  TsigFrame,
1609  ROOT.RooFit.ProjWData(
1610  ROOT.RooArgSet(
1611  self.qq),
1612  self.fitDataSigDeltaTPosfitDataSigDeltaTPos),
1613  ROOT.RooFit.LineWidth(3),
1614  ROOT.RooFit.LineColor(
1615  ROOT.kRed),
1616  ROOT.RooFit.Components("Tsigp"))
1617  TsigModel.plotOn(
1618  TsigFrame,
1619  ROOT.RooFit.ProjWData(
1620  ROOT.RooArgSet(
1621  self.qq),
1622  self.fitDataSigDeltaTPosfitDataSigDeltaTPos),
1623  ROOT.RooFit.LineWidth(3),
1624  ROOT.RooFit.LineColor(
1625  ROOT.kBlue + 1),
1626  ROOT.RooFit.LineStyle(7),
1627  ROOT.RooFit.Components("Tsigm"))
1628  TsigFrame.SetTitle("")
1629  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1630  TsigFrame.GetXaxis().SetTitleSize(0.05)
1631  TsigFrame.GetXaxis().SetLabelSize(0.045)
1632  TsigFrame.GetYaxis().SetTitleSize(0.05)
1633  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1634  TsigFrame.GetYaxis().SetLabelSize(0.045)
1635 
1636  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1637  c2.cd()
1638 
1639  self.TsigTsig.Print()
1640 
1641  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1642  Pad.SetLeftMargin(0.15)
1643  Pad.SetBottomMargin(0.15)
1644  Pad.Draw()
1645  Pad.cd()
1646  TsigFrame.Draw()
1647 
1648  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1649  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1650  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1651  curveMC.SetFillColor(ROOT.kWhite)
1652  curveMCB0.SetFillColor(ROOT.kWhite)
1653  curveMCAB0.SetFillColor(ROOT.kWhite)
1654 
1655  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1656  lMC.AddEntry(curveMC, "Both")
1657  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1658  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1659 
1660  lMC.SetTextSize(0.065)
1661  lMC.Draw()
1662 
1663  fitResSigPos.Print("v")
1664  c2.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMCPositiveDeltaT.pdf')
1665  c2.Clear()
1666 
1667  dotsMC = TsigFrame.getHist("h_fitDataSig")
1668  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1669  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1670 
1671  plotWithAsymmetry(
1672  TsigFrame,
1673  self.TsigPosNegTsigPosNeg,
1674  dotsMC,
1675  dotsMCB0,
1676  dotsMCAB0,
1677  curveMC,
1678  curveMCB0,
1679  curveMCAB0,
1680  "[ps]",
1681  "TsigFittedOnGenMCPositiveDeltaT",
1682  lMC,
1683  "",
1684  labelAsymmetry,
1685  True)
1686 
1687  # t_sig fit negative Delta t range
1688 
1689  # self.S.setConstant(ROOT.kTRUE)
1690  self.SS.setVal(0.)
1691  self.AA.setVal(0.)
1692  # self.A.setConstant(ROOT.kTRUE)
1693  self.AA.setConstant(ROOT.kFALSE)
1694 
1695  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg.Print()
1696 
1697  Tsigp = ROOT.RooGenericPdf(
1698  "Tsigp",
1699  "Tsigp",
1700  "(exp(-2*abs(Tsig)/Tau)/(Tau))*(1 + (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1701  ROOT.RooArgList(
1702  self.NormNorm,
1703  self.TsigTsig,
1704  self.TauTau,
1705  self.AA,
1706  self.SS,
1707  self.DMDM))
1708  Tsigm = ROOT.RooGenericPdf(
1709  "Tsigm",
1710  "Tsigm",
1711  "(exp(-2*abs(Tsig)/Tau)/(Tau))*(1 - (A - S*(DM*Tau))/(1+(Tau*DM)*(Tau*DM)))",
1712  ROOT.RooArgList(
1713  self.NormNorm,
1714  self.TsigTsig,
1715  self.TauTau,
1716  self.AA,
1717  self.SS,
1718  self.DMDM))
1719 
1720  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.qq)
1721  TsigModel.addPdf(Tsigp, "1")
1722  TsigModel.addPdf(Tsigm, "-1")
1723  TsigModel.Print()
1724 
1725  fitResSigNeg = TsigModel.fitTo(
1726  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg, ROOT.RooFit.Minos(
1727  ROOT.kFALSE), ROOT.RooFit.Extended(
1728  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1729  fitResSigNeg.Print("v")
1730 
1731  TsigFrame = self.TsigPosNegTsigPosNeg.frame()
1732  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg.plotOn(TsigFrame, ROOT.RooFit.Cut("Tsig < " + str(TsigFrame.GetXaxis().GetXmax())))
1733  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg.plotOn(TsigFrame,
1734  ROOT.RooFit.Cut("q==q::1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1735  ROOT.RooFit.MarkerColor(ROOT.kRed),
1736  ROOT.RooFit.Name("data_histB0sig"))
1737  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg.plotOn(TsigFrame,
1738  ROOT.RooFit.Cut("q==1 && Tsig < " + str(TsigFrame.GetXaxis().GetXmax())),
1739  ROOT.RooFit.MarkerColor(ROOT.kBlue + 1),
1740  ROOT.RooFit.Name("data_histAB0sig"))
1741 
1742  TsigModel.plotOn(TsigFrame, ROOT.RooFit.ProjWData(ROOT.RooArgSet(self.qq), self.fitDataSigDeltaTNegfitDataSigDeltaTNeg),
1743  ROOT.RooFit.LineWidth(3), ROOT.RooFit.LineColor(ROOT.kBlack), ROOT.RooFit.LineStyle(3))
1744  TsigModel.plotOn(
1745  TsigFrame,
1746  ROOT.RooFit.ProjWData(
1747  ROOT.RooArgSet(
1748  self.qq),
1749  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg),
1750  ROOT.RooFit.LineWidth(3),
1751  ROOT.RooFit.LineColor(
1752  ROOT.kRed),
1753  ROOT.RooFit.Components("Tsigp"))
1754  TsigModel.plotOn(
1755  TsigFrame,
1756  ROOT.RooFit.ProjWData(
1757  ROOT.RooArgSet(
1758  self.qq),
1759  self.fitDataSigDeltaTNegfitDataSigDeltaTNeg),
1760  ROOT.RooFit.LineWidth(3),
1761  ROOT.RooFit.LineColor(
1762  ROOT.kBlue + 1),
1763  ROOT.RooFit.LineStyle(7),
1764  ROOT.RooFit.Components("Tsigm"))
1765  TsigFrame.SetTitle("")
1766  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1767  TsigFrame.GetXaxis().SetTitleSize(0.05)
1768  TsigFrame.GetXaxis().SetLabelSize(0.045)
1769  TsigFrame.GetYaxis().SetTitleSize(0.05)
1770  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1771  TsigFrame.GetYaxis().SetLabelSize(0.045)
1772 
1773  c2 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1774  c2.cd()
1775 
1776  self.TsigTsig.Print()
1777 
1778  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1779  Pad.SetLeftMargin(0.15)
1780  Pad.SetBottomMargin(0.15)
1781  Pad.Draw()
1782  Pad.cd()
1783  TsigFrame.Draw()
1784  TsigFrame.Print()
1785 
1786  # lMC = ROOT.TLegend(0.7, 0.63, 0.9, 0.89)
1787 
1788  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1789  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1790  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1791  curveMC.SetFillColor(ROOT.kWhite)
1792  curveMCB0.SetFillColor(ROOT.kWhite)
1793  curveMCAB0.SetFillColor(ROOT.kWhite)
1794 
1795  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1796  lMC.AddEntry(curveMC, "Both")
1797  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1798  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1799 
1800  lMC.SetTextSize(0.065)
1801  lMC.Draw()
1802 
1803  fitResSigNeg.Print("v")
1804  c2.SaveAs("./" + 'B0_JPsiKs_TsigFittedOnGenMCNegativeDeltaT.pdf')
1805  c2.Clear()
1806 
1807  dotsMC = TsigFrame.getHist("h_fitDataSig")
1808  dotsMCB0 = TsigFrame.getHist("data_histB0sig")
1809  dotsMCAB0 = TsigFrame.getHist("data_histAB0sig")
1810 
1811  plotWithAsymmetry(
1812  TsigFrame,
1813  self.TsigPosNegTsigPosNeg,
1814  dotsMC,
1815  dotsMCB0,
1816  dotsMCAB0,
1817  curveMC,
1818  curveMCB0,
1819  curveMCAB0,
1820  "[ps]",
1821  "TsigFittedOnGenMCNegativeDeltaT",
1822  lMC,
1823  "",
1824  labelAsymmetry,
1825  True)
1826 
1827  # ----- Signal time fit qsig -------
1828 
1829  self.TauTau.setConstant(ROOT.kFALSE)
1830 
1831  self.fitDataSigfitDataSig.Print()
1832  Tsigp = ROOT.RooGenericPdf("Tsigp", "Tsigp", "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))",
1833  ROOT.RooArgList(self.NormNorm, self.TsigTsig, self.TauTau, self.DMDM))
1834  Tsigm = ROOT.RooGenericPdf("Tsigm", "Tsigm", "(exp(-1*abs(Tsig)/Tau)/(Norm*Tau))",
1835  ROOT.RooArgList(self.NormNorm, self.TsigTsig, self.TauTau, self.DMDM))
1836  TsigModel = ROOT.RooSimultaneous("TsigModel", "TsigModel", self.qsigqsig)
1837  TsigModel.addPdf(Tsigp, "1")
1838  TsigModel.addPdf(Tsigm, "-1")
1839  TsigModel.Print()
1840 
1841  fitResSig = TsigModel.fitTo(
1842  self.fitDataSigfitDataSig, ROOT.RooFit.Minos(
1843  ROOT.kFALSE), ROOT.RooFit.Extended(
1844  ROOT.kFALSE), ROOT.RooFit.NumCPU(1), ROOT.RooFit.Save())
1845  fitResSig.Print("v")
1846 
1847  TsigFrame = self.TsigTsig.frame()
1848  self.fitDataSigfitDataSig.plotOn(TsigFrame)
1849  self.fitDataSigfitDataSig.plotOn(
1850  TsigFrame,
1851  ROOT.RooFit.Cut("qsig==qsig::1"),
1852  ROOT.RooFit.MarkerColor(
1853  ROOT.kRed),
1854  ROOT.RooFit.Name("data_histB0sig"))
1855  self.fitDataSigfitDataSig.plotOn(
1856  TsigFrame,
1857  ROOT.RooFit.Cut("qsig==1"),
1858  ROOT.RooFit.MarkerColor(
1859  ROOT.kBlue + 1),
1860  ROOT.RooFit.Name("data_histAB0sig"))
1861 
1862  TsigModel.plotOn(
1863  TsigFrame,
1864  ROOT.RooFit.ProjWData(
1865  ROOT.RooArgSet(
1866  self.qsigqsig),
1867  self.fitDataSigfitDataSig),
1868  ROOT.RooFit.LineWidth(3),
1869  ROOT.RooFit.LineColor(
1870  ROOT.kBlack),
1871  ROOT.RooFit.LineStyle(3))
1872  TsigModel.plotOn(
1873  TsigFrame,
1874  ROOT.RooFit.ProjWData(
1875  ROOT.RooArgSet(
1876  self.qsigqsig),
1877  self.fitDataSigfitDataSig),
1878  ROOT.RooFit.LineWidth(3),
1879  ROOT.RooFit.LineColor(
1880  ROOT.kRed),
1881  ROOT.RooFit.Components("Tsigp"))
1882  TsigModel.plotOn(
1883  TsigFrame,
1884  ROOT.RooFit.ProjWData(
1885  ROOT.RooArgSet(
1886  self.qsigqsig),
1887  self.fitDataSigfitDataSig),
1888  ROOT.RooFit.LineWidth(3),
1889  ROOT.RooFit.LineColor(
1890  ROOT.kBlue + 1),
1891  ROOT.RooFit.LineStyle(7),
1892  ROOT.RooFit.Components("Tsigm"))
1893  TsigFrame.SetTitle("")
1894  TsigFrame.GetXaxis().SetTitle("#it{t}_{sig}^{gen} [ps]")
1895  TsigFrame.GetXaxis().SetTitleSize(0.05)
1896  TsigFrame.GetXaxis().SetLabelSize(0.045)
1897  TsigFrame.GetYaxis().SetTitleSize(0.05)
1898  TsigFrame.GetYaxis().SetTitleOffset(1.5)
1899  TsigFrame.GetYaxis().SetLabelSize(0.045)
1900 
1901  c3 = ROOT.TCanvas("c2", "c2", 1400, 1100)
1902  c3.cd()
1903 
1904  self.TsigTsig.Print()
1905 
1906  Pad = ROOT.TPad("p1", "p1", 0, 0, 1, 1, 0, 0, 0)
1907  Pad.SetLeftMargin(0.15)
1908  Pad.SetBottomMargin(0.15)
1909  Pad.Draw()
1910  Pad.cd()
1911  TsigFrame.Draw()
1912 
1913  curveMC = TsigFrame.getCurve("TsigModel_Norm[Tsig]")
1914  curveMCB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigp]")
1915  curveMCAB0 = TsigFrame.getCurve("TsigModel_Norm[Tsig]_Comp[Tsigm]")
1916  curveMC.SetFillColor(ROOT.kWhite)
1917  curveMCB0.SetFillColor(ROOT.kWhite)
1918  curveMCAB0.SetFillColor(ROOT.kWhite)
1919 
1920  lMC = ROOT.TLegend(0.73, 0.60, 0.89, 0.89)
1921  lMC.AddEntry(curveMC, "Both")
1922  lMC.AddEntry(curveMCAB0, '#bar{#it{B}}^{0}')
1923  lMC.AddEntry(curveMCB0, '#it{B}^{0}')
1924 
1925  lMC.SetTextSize(0.065)
1926  lMC.Draw()
1927 
1928  fitResSig.Print("v")
1929  c3.SaveAs("./" + 'B0_JPsiKs_TsigQsigFittedOnGenMC.pdf')
1930  c3.Clear()
1931 
1932  B0tagPercentage = self.B0sInTagSideB0sInTagSide / (self.B0barsInTagSideB0barsInTagSide + self.B0sInTagSideB0sInTagSide) * 100
1933  B0tagbarsPercentage = self.B0barsInTagSideB0barsInTagSide / (self.B0barsInTagSideB0barsInTagSide + self.B0sInTagSideB0sInTagSide) * 100
1934  Diff = B0tagbarsPercentage - B0tagPercentage
1935  print("B0 Tag Flavor: ")
1936  print(
1937  "B0s = ",
1938  f" {B0tagPercentage:0.4f}",
1939  "% B0bars = ",
1940  f" {B0tagbarsPercentage:0.4f}",
1941  "% Diff = ",
1942  f" {Diff:0.2f}",
1943  "%")
1944  B0signalPercentage = self.B0sInSignalSideB0sInSignalSide / (self.B0barsInSignalSideB0barsInSignalSide + self.B0sInSignalSideB0sInSignalSide) * 100
1945  B0signalbarsPercentage = self.B0barsInSignalSideB0barsInSignalSide / (self.B0barsInSignalSideB0barsInSignalSide + self.B0sInSignalSideB0sInSignalSide) * 100
1946  Diff = B0signalbarsPercentage - B0signalPercentage
1947  print("B0 Signal Flavor: ")
1948  print(
1949  "B0s = ",
1950  f" {B0signalPercentage:0.4f}",
1951  "% B0bars = ",
1952  f" {B0signalbarsPercentage:0.4f}",
1953  "% Diff = ",
1954  f" {Diff:0.2f}",
1955  "%")
1956  print("Total amount of events = ", self.nentriesnentries)
1957 
1958 
1959 fitDeltaTModule = fitDeltaT()
1960 validation_path.add_module(fitDeltaTModule)
1961 
1962 b2.process(validation_path)
1963 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