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