Belle II Software  release-08-01-10
CoGCalibration_utils.py
1 # !/usr/bin/env python
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 #
13 # evaluates the CoG corrections, create a localDB
14 # with the corrections and a root file to check
15 # the corrections
16 #
17 # usage: basf2 SVDCoGTimeCalibratinWithErrorImporter localDB filename
18 # localDB = name of the local DB folder
19 # filename = single root file, or file with the list of reconstructed files
20 #
21 # this script can be launched with launch_calibration_cog.sh in the
22 # B2SVD project, svd_CoGHitTime_calibration repository
23 #
24 
25 
26 import basf2 as b2
27 from ROOT import Belle2, TFile, TFitResultPtr, TH1F, TH2D, TH2F
28 from ROOT import TF1, gDirectory, gROOT
29 import math
30 from ROOT.Belle2 import SVDCoGCalibrationFunction
31 
32 
33 svd_recoDigits = "SVDRecoDigitsFromTracks"
34 cdc_Time0 = "EventT0"
35 svd_Clusters = "SVDClustersFromTracks"
36 
37 gROOT.SetBatch(True)
38 
39 # mode = True
40 
41 
43  """
44  Python class used for evaluating the CoG corrections, create a localDB,
45  creating a localDB with the corrections and a root file to check the corrections
46  """
47 
49  """
50  Function that allows to set if apply the CDC latency correction
51 
52  parameters:
53  mode (bool):
54  - if True -> not apply correction
55  - if False -> apply correction
56  """
57 
58  self.notApplyCDCLatencyCorrectionnotApplyCDCLatencyCorrection = mode
59  print("Not Correct for CDC latency: " + str(mode) + " " + str(self.notApplyCDCLatencyCorrectionnotApplyCDCLatencyCorrection))
60 
61  def fillLists(self, mode_byte_object, svdClusters_rel_RecoTracks_cl):
62  """
63  Function that fill the lists needed for the CoG corrections
64 
65  parameters:
66  mode_byte_object (modeByte): modeByte that contains the information about the TB
67  svdClusters_rel_RecoTracks_cl (SVDCluster): cluster related to tracks
68  """
69 
70  timeCluster = svdClusters_rel_RecoTracks_cl.getClsTime()
71  snrCluster = svdClusters_rel_RecoTracks_cl.getSNR()
72  layerCluster = svdClusters_rel_RecoTracks_cl.getSensorID().getLayerNumber()
73  layerIndex = layerCluster - 3
74  sensorCluster = svdClusters_rel_RecoTracks_cl.getSensorID().getSensorNumber()
75  sensorIndex = sensorCluster - 1
76  ladderCluster = svdClusters_rel_RecoTracks_cl.getSensorID().getLadderNumber()
77  ladderIndex = ladderCluster - 1
78  sideCluster = svdClusters_rel_RecoTracks_cl.isUCluster()
79  if sideCluster:
80  sideIndex = 1
81  else:
82  sideIndex = 0
83 
84  hasTimezero = self.cdcEventT0cdcEventT0.hasEventT0()
85  # print("Time: " + str(hasTimezero))
86  if hasTimezero:
87  TBClusters = mode_byte_object.getTriggerBin()
88  TBIndex = ord(TBClusters)
89  tZero = self.cdcEventT0cdcEventT0.getEventT0()
90  # tZero_err = self.cdcEventT0.getEventT0Uncertainty()
91  # tZero_err = 5.1
92  tZeroSync = tZero - 4000./509 * (3 - TBIndex)
93  et0 = self.EventT0HistEventT0Hist
94  et0.Fill(tZeroSync)
95  # print(str(tZero_err))
96 
97  resHist = self.resListresList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex]
98  resHist.Fill(timeCluster - tZeroSync)
99  spHist = self.spListspList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex]
100  # spHist.Fill(timeCluster, 1.3*timeCluster - 50 + random.gauss(0,10))
101  spHist.Fill(timeCluster, tZeroSync)
102  cogHist = self.cogListcogList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex]
103  cogHist.Fill(timeCluster)
104  cdcHist = self.cdcListcdcList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex]
105  cdcHist.Fill(tZeroSync)
106  snrHist = self.snrListsnrList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex]
107  snrHist.Fill(snrCluster)
108 
109  self.nListnList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex] += 1
110  self.sumCOGListsumCOGList[layerIndex][ladderIndex][sensorIndex][sideIndex][TBIndex] += timeCluster
111 
112 
113  self.NTOTNTOT = self.NTOTNTOT + 1
114 
115  def set_localdb(self, localDB):
116  """
117  Function that allows to set the localDB
118 
119  parameters:
120  localDB (str): Name of the localDB used
121  """
122 
123  self.localdblocaldb = localDB
124 
125  def initialize(self):
126  """
127  Initialize object (histograms, lists, ...) used by the class
128  """
129 
130 
131  self.outputFileNameoutputFileName = "CoGCorrectionMonitor_" + self.localdblocaldb + ".root"
132 
133 
135  self.resListresList = []
136 
137  self.spListspList = []
138 
139  self.cogListcogList = []
140 
141  self.cdcListcdcList = []
142 
143  self.snrListsnrList = []
144 
145  self.nListnList = []
146 
147  self.sumCOGListsumCOGList = []
148 
150 
151  self.EvtEvt = 0
152  for layer in geoCache.getLayers(Belle2.VXD.SensorInfoBase.SVD):
153  layerList0 = []
154  layerList1 = []
155  layerList2 = []
156  layerList3 = []
157  layerList4 = []
158  layerList5 = []
159  layerList6 = []
160  layerList7 = []
161  layerList8 = []
162  self.resListresList.append(layerList0)
163  self.spListspList.append(layerList1)
164  self.cogListcogList.append(layerList2)
165  self.cdcListcdcList.append(layerList3)
166  self.snrListsnrList.append(layerList4)
167  self.nListnList.append(layerList8)
168  self.sumCOGListsumCOGList.append(layerList7)
169  # layerNumber = layer.getLayerNumber()
170  for ladder in geoCache.getLadders(layer):
171  ladderList0 = []
172  ladderList1 = []
173  ladderList2 = []
174  ladderList3 = []
175  ladderList4 = []
176  ladderList5 = []
177  ladderList6 = []
178  ladderList7 = []
179  ladderList8 = []
180  layerList0.append(ladderList0)
181  layerList1.append(ladderList1)
182  layerList2.append(ladderList2)
183  layerList3.append(ladderList3)
184  layerList4.append(ladderList4)
185  layerList5.append(ladderList5)
186  layerList6.append(ladderList6)
187  layerList7.append(ladderList7)
188  layerList8.append(ladderList8)
189  # ladderNumber = ladder.getLadderNumber()
190  for sensor in geoCache.getSensors(ladder):
191  sensorList0 = []
192  sensorList1 = []
193  sensorList2 = []
194  sensorList3 = []
195  sensorList4 = []
196  sensorList5 = []
197  sensorList6 = []
198  sensorList7 = []
199  sensorList8 = []
200  ladderList0.append(sensorList0)
201  ladderList1.append(sensorList1)
202  ladderList2.append(sensorList2)
203  ladderList3.append(sensorList3)
204  ladderList4.append(sensorList4)
205  ladderList5.append(sensorList5)
206  ladderList6.append(sensorList6)
207  ladderList7.append(sensorList7)
208  ladderList8.append(sensorList8)
209  # sensorNumber = sensor.getSensorNumber()
210  for side in range(2):
211  sideList0 = []
212  sideList1 = []
213  sideList2 = []
214  sideList3 = []
215  sideList4 = []
216  sideList5 = []
217  sideList6 = []
218  sideList7 = []
219  sideList8 = []
220  sensorList0.append(sideList0)
221  sensorList1.append(sideList1)
222  sensorList2.append(sideList2)
223  sensorList3.append(sideList3)
224  sensorList4.append(sideList4)
225  sensorList5.append(sideList5)
226  sensorList6.append(sideList6)
227  sensorList7.append(sideList7)
228  sensorList8.append(sideList8)
229 
230  for i in geoCache.getLayers(Belle2.VXD.SensorInfoBase.SVD):
231  layerN = i.getLayerNumber()
232  li = layerN - 3
233  for j in geoCache.getLadders(i):
234  ladderN = j.getLadderNumber()
235  ldi = ladderN - 1
236  for k in geoCache.getSensors(j):
237  sensorN = k.getSensorNumber()
238  si = sensorN - 1
239  for s in range(2):
240  for t in range(4):
241  self.resListresList[li][ldi][si][s].append(
242  TH1F("res" + "_" + str(k) + "." + str(s) + "." + str(t), " ", 200, -100, 100))
243  self.spListspList[li][ldi][si][s].append(
244  TH2D("sp" + "_" + str(k) + "." + str(s) + "." + str(t), " ", 300, -150, 150, 300, -150, 150))
245  self.cogListcogList[li][ldi][si][s].append(
246  TH1F("cog" + "_" + str(k) + "." + str(s) + "." + str(t), " ", 200, -100, 100))
247  self.cdcListcdcList[li][ldi][si][s].append(
248  TH1F("cdc" + "_" + str(k) + "." + str(s) + "." + str(t), " ", 200, -100, 100))
249  self.snrListsnrList[li][ldi][si][s].append(
250  TH1F("snr" + "_" + str(k) + "." + str(s) + "." + str(t), " ", 100, 0, 100))
251  self.nListnList[li][ldi][si][s].append(0)
252  self.sumCOGListsumCOGList[li][ldi][si][s].append(0)
253 
254  self.EventT0HistEventT0Hist = TH1F("EventT0", " ", 200, -100, 100)
255 
256  self.AlphaUTBAlphaUTB = TH2F("alphaVsTB_U", " ", 400, 0.5, 2, 4, 0, 4)
257  self.AlphaUTBAlphaUTB.GetXaxis().SetTitle("alpha")
258  self.AlphaUTBAlphaUTB.GetYaxis().SetTitle("trigger bin")
259 
260  self.AlphaVTBAlphaVTB = TH2F("alphaVsTB_V", " ", 400, 0.5, 2, 4, 0, 4)
261  self.AlphaVTBAlphaVTB.GetXaxis().SetTitle("alpha")
262  self.AlphaVTBAlphaVTB.GetYaxis().SetTitle("trigger bin")
263 
264  self.BetaUTBBetaUTB = TH2F("betaVsTB_U", " ", 200, -100, 100, 4, 0, 4)
265  self.BetaUTBBetaUTB.GetXaxis().SetTitle("beta (ns)")
266  self.BetaUTBBetaUTB.GetYaxis().SetTitle("trigger bin")
267 
268  self.BetaVTBBetaVTB = TH2F("betaVsTB_V", " ", 200, -100, 100, 4, 0, 4)
269  self.BetaVTBBetaVTB.GetXaxis().SetTitle("beta (ns)")
270  self.BetaVTBBetaVTB.GetYaxis().SetTitle("trigger bin")
271 
272 
273  self.MeanHistVTBMeanHistVTB = TH2F("meanHistVsTB_V", " ", 100, -10, 10, 4, 0, 4)
274  self.MeanHistVTBMeanHistVTB.GetXaxis().SetTitle("distribution mean (ns)")
275  self.MeanHistVTBMeanHistVTB.GetYaxis().SetTitle("trigger bin")
276 
277  self.MeanHistUTBMeanHistUTB = TH2F("meanHistVsTB_U", " ", 100, -10, 10, 4, 0, 4)
278  self.MeanHistUTBMeanHistUTB.GetXaxis().SetTitle("distribution mean (ns)")
279  self.MeanHistUTBMeanHistUTB.GetYaxis().SetTitle("trigger bin")
280 
281  self.RMSHistVTBRMSHistVTB = TH2F("rmsHistVsTB_V", " ", 100, 0, 10, 4, 0, 4)
282  self.RMSHistVTBRMSHistVTB.GetXaxis().SetTitle("distribution RMS (ns)")
283  self.RMSHistVTBRMSHistVTB.GetYaxis().SetTitle("trigger bin")
284 
285  self.RMSHistUTBRMSHistUTB = TH2F("rmsHistVsTB_U", " ", 100, 0, 10, 4, 0, 4)
286  self.RMSHistUTBRMSHistUTB.GetXaxis().SetTitle("distribution RMS (ns)")
287  self.RMSHistUTBRMSHistUTB.GetYaxis().SetTitle("trigger bin")
288 
289  self.MeanFitVTBMeanFitVTB = TH2F("meanFitVsTB_V", " ", 100, -10, 10, 4, 0, 4)
290  self.MeanFitVTBMeanFitVTB.GetXaxis().SetTitle("fit mean (ns)")
291  self.MeanFitVTBMeanFitVTB.GetYaxis().SetTitle("trigger bin")
292 
293  self.MeanFitUTBMeanFitUTB = TH2F("meanFitVsTB_U", " ", 100, -10, 10, 4, 0, 4)
294  self.MeanFitUTBMeanFitUTB.GetXaxis().SetTitle("fit mean (ns)")
295  self.MeanFitUTBMeanFitUTB.GetYaxis().SetTitle("trigger bin")
296 
297  self.RMSFitUTBRMSFitUTB = TH2F("rmsFitVsTB_U", " ", 100, 0, 10, 4, 0, 4)
298  self.RMSFitUTBRMSFitUTB.GetXaxis().SetTitle("fit sigma (ns)")
299  self.RMSFitUTBRMSFitUTB.GetYaxis().SetTitle("trigger bin")
300 
301  self.RMSFitVTBRMSFitVTB = TH2F("rmsFitVsTB_V", " ", 100, 0, 10, 4, 0, 4)
302  self.RMSFitVTBRMSFitVTB.GetXaxis().SetTitle("fit sigma (ns)")
303  self.RMSFitVTBRMSFitVTB.GetYaxis().SetTitle("trigger bin")
304 
305 
306  self.gausgaus = TF1("gaus", 'gaus(0)', -150, 100)
307 
308  self.NTOTNTOT = 0
309 
310  def event(self):
311  """
312  Function that allows to cicle on the events
313  """
314  svd_evt_info = Belle2.PyStoreObj('SVDEventInfo')
315  mode_byte = svd_evt_info.getModeByte()
316  # timeClusterU = 0
317  # timeClusterV = 0
318  # sideIndex = 0
319  # TBIndexU = 0
320  # TBIndexV = 0
321 
322  self.EvtEvt = self.EvtEvt + 1
323 
324 
325  self.cdcEventT0cdcEventT0 = Belle2.PyStoreObj(cdc_Time0)
326  svdCluster_list = Belle2.PyStoreArray(svd_Clusters)
327  # svdRecoDigit_list = Belle2.PyStoreArray(svd_recoDigits)
328 
329  for svdCluster in svdCluster_list:
330  # svdRecoDigit = svdCluster.getRelatedTo(svd_recoDigits)
331  self.fillListsfillLists(mode_byte, svdCluster)
332 
333  def terminate(self):
334  """
335  Terminates te class and produces the output rootfile
336  """
337 
338  tfile = TFile(self.outputFileNameoutputFileName, 'recreate')
340 
342 
343  timeCal = SVDCoGCalibrationFunction()
344  # Bias and Scale
345  tbBias = [-50, -50, -50, -50]
346  tbScale = [1, 1, 1, 1]
347  tbBias_err = [1, 1, 1, 1]
348  tbScale_err = [1, 1, 1, 1]
349  # tbCovScaleBias = [1, 1, 1, 1]
350 
351  TCOGMEAN = 0
352  T0MEAN = 0
353 
354  '''
355  geoCache = Belle2.VXD.GeoCache.getInstance()
356  for layer in geoCache.getLayers(Belle2.VXD.SensorInfoBase.SVD):
357  layerNumber = layer.getLayerNumber()
358  li = layerNumber - 3
359  for ladder in geoCache.getLadders(layer):
360  ladderNumber = ladder.getLadderNumber()
361  ldi = ladderNumber - 1
362  for sensor in geoCache.getSensors(ladder):
363  sensorNumber = sensor.getSensorNumber()
364  si = sensorNumber - 1
365  for side in range(2):
366  for tb in range(4):
367  n = self.nList[li][ldi][si][side][tb]
368  NTOT += n
369  '''
371  gDirectory.mkdir("plots")
372  gDirectory.cd("plots")
373  for layer in geoCache.getLayers(Belle2.VXD.SensorInfoBase.SVD):
374  layerNumber = layer.getLayerNumber()
375  li = layerNumber - 3
376  gDirectory.mkdir("layer" + str(layer))
377  gDirectory.cd("layer" + str(layer))
378  for ladder in geoCache.getLadders(layer):
379  ladderNumber = ladder.getLadderNumber()
380  ldi = ladderNumber - 1
381  for sensor in geoCache.getSensors(ladder):
382  sensorNumber = sensor.getSensorNumber()
383  si = sensorNumber - 1
384  for side in range(2):
385  for tb in range(4):
386  # Resolution distribution Histograms with Gaussian Fit
387  res = self.resListresList[li][ldi][si][side][tb]
388  fitResult = int(TFitResultPtr(res.Fit(self.gausgaus, "R")))
389 
390  if res.GetEntries() > 5:
391  if side == 1:
392  self.MeanHistUTBMeanHistUTB.Fill(res.GetMean(), tb)
393  self.RMSHistUTBRMSHistUTB.Fill(res.GetRMS(), tb)
394  if fitResult > -1:
395  self.MeanFitUTBMeanFitUTB.Fill(self.gausgaus.GetParameter(1), tb)
396  self.RMSFitUTBRMSFitUTB.Fill(self.gausgaus.GetParameter(2), tb)
397  else:
398  self.MeanHistVTBMeanHistVTB.Fill(res.GetMean(), tb)
399  self.RMSHistVTBRMSHistVTB.Fill(res.GetRMS(), tb)
400  if fitResult > -1:
401  self.MeanFitVTBMeanFitVTB.Fill(self.gausgaus.GetParameter(1), tb)
402  self.RMSFitVTBRMSFitVTB.Fill(self.gausgaus.GetParameter(2), tb)
403 
404  res.Write()
405  # COG Distribution Histograms
406  cog = self.cogListcogList[li][ldi][si][side][tb]
407  cog.Write()
408  # CDC EventT0 Distribution Histograms
409  cdc = self.cdcListcdcList[li][ldi][si][side][tb]
410  cdc.Write()
411  # SNR Distribution Histograms
412  snr = self.snrListsnrList[li][ldi][si][side][tb]
413  # snrMean = snr.GetMean()
414  snr.Write()
415  # ScatterPlot Histograms with Linear Fit
416  sp = self.spListspList[li][ldi][si][side][tb]
417  # covscalebias = sp.GetCovariance()
418  pfxsp = sp.ProfileX()
419  sp.Write()
420  pfxsp.Write()
421 
422  if sp.GetRMS() != 0:
423  m = sp.GetCovariance() / pow(sp.GetRMS(1), 2)
424  # m = sp.GetCovariance()/cog.GetRMS()
425  m_err = 2 / pow(sp.GetRMS(), 3) * sp.GetRMSError() * sp.GetCovariance()
426  q = sp.GetMean(2) - m * sp.GetMean(1)
427  q_err = math.sqrt(pow(sp.GetMeanError(2), 2) +
428  pow(m * sp.GetMeanError(1), 2) + pow(m_err * sp.GetMean(1), 2))
429  else:
430  m = 1
431  m_err = 0
432  q = 0
433  q_err = 0
434 
435  if side == 1:
436  self.AlphaUTBAlphaUTB.Fill(m, tb)
437  self.BetaUTBBetaUTB.Fill(q, tb)
438  else:
439  self.AlphaVTBAlphaVTB.Fill(m, tb)
440  self.BetaVTBBetaVTB.Fill(q, tb)
441 
442  n = self.nListnList[li][ldi][si][side][tb]
443 
444  tbBias[tb] = q
445  tbScale[tb] = m
446  tbBias_err[tb] = q_err
447  tbScale_err[tb] = m_err
448  TCOGMEAN += n * (m * self.sumCOGListsumCOGList[li][ldi][si][side][tb] / n + q) / self.NTOTNTOT
449 
450  T0MEAN = self.EventT0HistEventT0Hist.GetMean()
451  '''
452  print(
453  "Mean of the CoG corrected distribution: " +
454  str(TCOGMEAN) +
455  " Mean of the T0 distribution: " +
456  str(T0MEAN))
457  '''
458  if not self.notApplyCDCLatencyCorrectionnotApplyCDCLatencyCorrection:
459  tbBias[0] = tbBias[0] - T0MEAN
460  tbBias[1] = tbBias[1] - T0MEAN
461  tbBias[2] = tbBias[2] - T0MEAN
462  tbBias[3] = tbBias[3] - T0MEAN
463 
464  timeCal.set_bias(tbBias[0], tbBias[1], tbBias[2], tbBias[3])
465  timeCal.set_scale(tbScale[0], tbScale[1], tbScale[2], tbScale[3])
466  print("setting CoG calibration for " + str(layerNumber) + "." + str(ladderNumber) + "." + str(sensorNumber))
467  payload.set(layerNumber, ladderNumber, sensorNumber, bool(side), 1, timeCal)
468  gDirectory.cd("../")
469 
470  gDirectory.cd("../")
471  self.AlphaUTBAlphaUTB.Write()
472  self.AlphaVTBAlphaVTB.Write()
473  self.BetaUTBBetaUTB.Write()
474  self.BetaVTBBetaVTB.Write()
475  self.MeanFitUTBMeanFitUTB.Write()
476  self.MeanFitVTBMeanFitVTB.Write()
477  self.RMSFitUTBRMSFitUTB.Write()
478  self.RMSFitVTBRMSFitVTB.Write()
479  self.MeanHistUTBMeanHistUTB.Write()
480  self.MeanHistVTBMeanHistVTB.Write()
481  self.RMSHistUTBRMSHistUTB.Write()
482  self.RMSHistVTBRMSHistVTB.Write()
483  self.EventT0HistEventT0Hist.Write()
484 
485  Belle2.Database.Instance().storeData(Belle2.SVDCoGTimeCalibrations.name, payload, iov)
486 
487  tfile.Close()
static IntervalOfValidity always()
Function that returns an interval of validity that is always valid, c.f.
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
base class for calibrations classes
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
MeanHistVTB
mean of the residuals distribution vs TB, for V side
notApplyCDCLatencyCorrection
parameter that allows to apply or not the CDC latency correction
RMSFitUTB
RMS of the residuals distribution vs TB, for U side (from gaussian fit)
resList
lists used to create the histograms for each TB : residuals
RMSFitVTB
RMS of the residuals distribution vs TB, for V side (from gaussian fit)
MeanFitUTB
mean of the residuals distribution vs TB, for U side (from gaussian fit)
def fillLists(self, mode_byte_object, svdClusters_rel_RecoTracks_cl)
MeanFitVTB
mean of the residuals distribution vs TB, for V side (from gaussian fit)
MeanHistUTB
mean of the residuals distribution vs TB, for U side
RMSHistVTB
RMS of the residuals distribution vs TB, for V side.
RMSHistUTB
RMS of the residuals distribution vs TB, for U side.
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42