Belle II Software  release-06-01-15
cdst_chi2ModuleT0calibration.py
1 #!/usr/bin/env python
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # --------------------------------------------------------------------------------
13 # Calibrate module T0 with Bhabha's or dimuons using time difference between slots
14 #
15 # usage: basf2 cdst_chi2ModuleT0calibration.py experiment runFirst runLast
16 # job: bsub -q l "basf2 cdst_chi2ModuleT0calibration.py experiment runFirst runLast"
17 #
18 # note: runLast is inclusive
19 # --------------------------------------------------------------------------------
20 
21 import basf2 as b2
22 from ROOT import Belle2
23 from ROOT import TFile, TH1F, TH2F, TF1, TMatrixDSym
24 from array import array
25 import math
26 import sys
27 import glob
28 import os
29 
30 # ----- those need to be adjusted before running --------------------------------------
31 #
32 sampleType = 'bhabha' # sample type: 'bhabha' or 'dimuon'
33 data_dir = '/group/belle2/dataprod/Data/release-03-02-02/DB00000635/proc00000009_nofilter'
34 skim_dir = 'skim/hlt_bhabha/cdst/sub00/'
35 globalTag = 'data_reprocessing_prompt' # base global tag
36 stagingTags = ['staging_data_reprocessing'] # list of tags with new calibration
37 localDB = [] # list of local databases with new calibration
38 minEntries = 10 # minimal number of histogram entries required to perform a fit
39 output_dir = 'moduleT0' # main output folder
40 #
41 # -------------------------------------------------------------------------------------
42 
43 # Argument parsing
44 argvs = sys.argv
45 if len(argvs) < 4:
46  print("usage: basf2", argvs[0], "experiment runFirst runLast")
47  sys.exit()
48 experiment = int(argvs[1])
49 run_first = int(argvs[2])
50 run_last = int(argvs[3])
51 
52 expNo = 'e' + '{:0=4d}'.format(experiment)
53 
54 # Make list of files
55 files = []
56 for run in range(run_first, run_last + 1):
57  runNo = 'r' + '{:0=5d}'.format(run)
58  for typ in ['4S', 'Continuum', 'Scan']:
59  folder = data_dir + '/' + expNo + '/' + typ + '/' + runNo + '/' + skim_dir
60  files += glob.glob(folder + '/cdst.*.root')
61 if len(files) == 0:
62  b2.B2ERROR('No cdst files found')
63  sys.exit()
64 
65 # Output folder
66 method = 'chi2'
67 output_folder = output_dir + '/' + expNo + '/' + sampleType + '/' + method
68 if not os.path.isdir(output_folder):
69  os.makedirs(output_folder)
70  print('New folder created: ' + output_folder)
71 
72 # Output file name
73 fileName = output_folder + '/moduleT0-' + expNo + '-'
74 run1 = 'r' + '{:0=5d}'.format(run_first)
75 run2 = 'r' + '{:0=5d}'.format(run_last)
76 fileName += run1 + '_to_' + run2 + '.root'
77 print('Output file:', fileName)
78 
79 
80 class Mask_BS13d(b2.Module):
81  ''' exclude (mask-out) BS 13d '''
82 
83  def event(self):
84  ''' event processing '''
85 
86  for digit in Belle2.PyStoreArray('TOPDigits'):
87  if digit.getModuleID() == 13 and digit.getBoardstackNumber() == 3:
88  digit.setHitQuality(Belle2.TOPDigit.c_Junk)
89 
90 
91 class ModuleT0cal(b2.Module):
92  ''' module T0 calibrator using chi2 minimization of time differences between slots '''
93 
94  def initialize(self):
95  ''' initialize and book histograms '''
96 
97  Belle2.PyStoreArray('TOPTimeZeros').isRequired()
98 
99 
100  self.db_moduleT0db_moduleT0 = Belle2.PyDBObj('TOPCalModuleT0')
101 
102 
103  self.h_deltaT0h_deltaT0 = {}
104 
105 
106  self.h_slotPairsh_slotPairs = TH2F("slots", "slot pairs: number of entries",
107  16, 0.5, 16.5, 16, 0.5, 16.5)
108  self.h_slotPairsh_slotPairs.SetXTitle("first slot number")
109  self.h_slotPairsh_slotPairs.SetYTitle("second slot number")
110 
111 
112  self.h_slotPairs_acch_slotPairs_acc = TH2F("slots_acc",
113  "slot pairs: #chi^{2} /NDF of " +
114  "successfully fitted time differences",
115  16, 0.5, 16.5, 16, 0.5, 16.5)
116  self.h_slotPairs_acch_slotPairs_acc.SetXTitle("first slot number")
117  self.h_slotPairs_acch_slotPairs_acc.SetYTitle("second slot number")
118 
119 
120  self.h_relModuleT0h_relModuleT0 = TH1F("relModuleT0", "Module T0 relative to calibration",
121  16, 0.5, 16.5)
122  self.h_relModuleT0h_relModuleT0.SetXTitle("slot number")
123  self.h_relModuleT0h_relModuleT0.SetYTitle("module T0 residual [ns]")
124 
125 
126  self.h_moduleT0h_moduleT0 = TH1F("moduleT0", "Module T0", 16, 0.5, 16.5)
127  self.h_moduleT0h_moduleT0.SetXTitle("slot number")
128  self.h_moduleT0h_moduleT0.SetYTitle("module T0 [ns]")
129 
130  def sortedTimeZeros(self, unsortedPyStoreArray):
131  ''' sorting timeZeros according to slot number '''
132 
133  # first convert to a python-list to be able to sort
134  py_list = [x for x in unsortedPyStoreArray]
135  # sort via a hierachy of sort keys
136  return sorted(
137  py_list,
138  key=lambda x: (
139  x.getModuleID(),
140  )
141  )
142 
143  def event(self):
144  ''' event processor: fill histograms for events with two entries in timeZeros '''
145 
146  timeZeros = Belle2.PyStoreArray('TOPTimeZeros')
147  if timeZeros.getEntries() != 2:
148  return
149 
150  tZeros = self.sortedTimeZerossortedTimeZeros(timeZeros)
151 
152  slot1 = tZeros[0].getModuleID()
153  slot2 = tZeros[1].getModuleID()
154  self.h_slotPairsh_slotPairs.Fill(slot1, slot2)
155  if slot1 not in self.h_deltaT0h_deltaT0:
156  self.h_deltaT0h_deltaT0[slot1] = {}
157  if slot2 not in self.h_deltaT0h_deltaT0[slot1]:
158  name = 'deltaT0_' + str(slot1) + '-' + str(slot2)
159  title = 'time difference: slot ' + str(slot1) + ' - slot ' + str(slot2)
160  self.h_deltaT0h_deltaT0[slot1][slot2] = TH1F(name, title, 200, -5, 5)
161  self.h_deltaT0h_deltaT0[slot1][slot2].SetXTitle('time difference [ns]')
162  dt = tZeros[0].getTime() - tZeros[1].getTime()
163  self.h_deltaT0h_deltaT0[slot1][slot2].Fill(dt)
164 
165  def singleFit(self, h):
166  ''' perform a fit of time difference distribution '''
167 
168  xmin = h.GetXaxis().GetXmin()
169  xmax = h.GetXaxis().GetXmax()
170  fun = TF1("fun", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]", xmin, xmax)
171  fun.SetParameter(0, h.GetMaximum())
172  fun.SetParameter(1, h.GetBinCenter(h.GetMaximumBin()))
173  fun.SetParameter(2, 0.2)
174  fun.SetParameter(3, 0)
175  return h.Fit(fun, "LERSQ")
176 
177  def fitHistograms(self, minEntries):
178  ''' fit all histograms and store the results into lists '''
179 
180 
181  self.delT0delT0 = []
182 
183  self.sigmasigma = []
184 
185  self.AA = []
186 
187  self.h_outh_out = []
188 
189  for slot1 in range(1, 10):
190  if slot1 not in self.h_deltaT0h_deltaT0:
191  continue
192  for slot2 in range(slot1 + 7, slot1 + 10):
193  if slot2 not in self.h_deltaT0h_deltaT0[slot1]:
194  continue
195  h = self.h_deltaT0h_deltaT0[slot1][slot2]
196  self.h_outh_out.append(h)
197  if h.GetEntries() > minEntries:
198  result = self.singleFitsingleFit(h)
199  if int(result) == 0:
200  chi2 = result.Chi2() / result.Ndf()
201  self.h_slotPairs_acch_slotPairs_acc.SetBinContent(slot1, slot2, chi2)
202  self.delT0delT0.append(result.Parameter(1))
203  self.sigmasigma.append(result.ParError(1))
204  a = [0 for i in range(16)]
205  a[slot1 - 1] = 1
206  a[slot2 - 1] = -1
207  self.AA.append(a)
208 
209  def minimize(self):
210  '''
211  Minimization procedure. For the method see NIM A 586 (2008) 174-179, sect. 2.2.
212  '''
213 
214  # append the bound (sum of all calibration constants equals to 0)
215  self.delT0delT0.append(0.0)
216  self.sigmasigma.append(0.001) # arbitrary but small compared to calibration precision
217  a = [1 for i in range(16)]
218  self.AA.append(a)
219  m = len(self.delT0delT0)
220 
221  # for i, a in enumerate(self.A):
222  # print(a, round(self.delT0[i], 3), round(self.sigma[i], 3))
223 
224  # construct the matrix of a linear system of equations
225  B = TMatrixDSym(16)
226  for i in range(16):
227  for j in range(16):
228  for k in range(m):
229  B[i][j] += self.AA[k][i] * self.AA[k][j] / self.sigmasigma[k]**2
230 
231  # invert the matrix
232  det = array('d', [0])
233  B.Invert(det)
234  if det[0] == 0:
235  b2.B2ERROR("Matrix inversion failed")
236  return False
237 
238  # construct the right side of a linear system of equations
239  b = [0.0 for i in range(16)]
240  for i in range(16):
241  for k in range(m):
242  b[i] += self.AA[k][i] * self.delT0delT0[k] / self.sigmasigma[k]**2
243 
244  # solve for unknown relative module T0
245  x = [0 for i in range(16)] # module T0's
246  e = [0 for i in range(16)] # uncertainties on module T0
247  for i in range(16):
248  for j in range(16):
249  x[i] += B[i][j] * b[j]
250  e[i] = math.sqrt(B[i][i])
251 
252  # calculate chi^2 and ndf
253  chi2 = 0
254  for k in range(m):
255  s = 0
256  for i in range(16):
257  s += self.AA[k][i] * x[i]
258  chi2 += ((s - self.delT0delT0[k]) / self.sigmasigma[k])**2
259  ndf = m - 16
260  chi_ndf = 'chi^2/NDF = ' + str(round(chi2, 2)) + '/' + str(ndf)
261  self.h_relModuleT0h_relModuleT0.SetTitle(self.h_relModuleT0h_relModuleT0.GetTitle() + ' (' + chi_ndf + ')')
262 
263  # store relative module T0 in a histogram
264  for i in range(16):
265  self.h_relModuleT0h_relModuleT0.SetBinContent(i+1, x[i])
266  self.h_relModuleT0h_relModuleT0.SetBinError(i+1, e[i])
267 
268  # calculate arithmetic average of current calibration constants
269  T0 = [0.0 for i in range(16)]
270  average = 0.0
271  num = 0
272  for i in range(16):
273  if self.db_moduleT0db_moduleT0.isCalibrated(i+1):
274  T0[i] = self.db_moduleT0db_moduleT0.getT0(i+1)
275  average += T0[i]
276  num += 1
277  if num > 0:
278  average /= num
279 
280  # add current calibration to relative one and subtract the average, then store
281  for i in range(16):
282  self.h_moduleT0h_moduleT0.SetBinContent(i+1, x[i] + T0[i] - average)
283  self.h_moduleT0h_moduleT0.SetBinError(i+1, e[i])
284 
285  b2.B2RESULT("Module T0 calibration constants successfully determined, " + chi_ndf)
286  return True
287 
288  def terminate(self):
289  ''' fit histograms, minimize and write results to file '''
290 
291  # fit time difference distributions
292  self.fitHistogramsfitHistograms(minEntries)
293 
294  # minimize
295  OK = self.minimizeminimize()
296 
297  # open output file
298  file = TFile.Open(fileName, 'recreate')
299 
300  # write histograms and close the file
301  self.h_slotPairsh_slotPairs.Write()
302  self.h_slotPairs_acch_slotPairs_acc.Write()
303  if OK:
304  self.h_relModuleT0h_relModuleT0.Write()
305  self.h_moduleT0h_moduleT0.Write()
306  for h in self.h_outh_out:
307  h.Write()
308  file.Close()
309 
310  b2.B2RESULT("Output written to " + fileName)
311 
312 
313 # Database
314 b2.use_central_database(globalTag)
315 for tag in stagingTags:
316  b2.use_central_database(tag)
317 for db in localDB:
318  if os.path.isfile(db):
319  b2.use_local_database(db, invertLogging=True)
320  else:
321  b2.B2ERROR(db + ": local database not found")
322  sys.exit()
323 
324 # Create paths
325 main = b2.create_path()
326 
327 # Input (cdst files)
328 main.add_module('RootInput', inputFileNames=files)
329 
330 # Initialize TOP geometry parameters
331 main.add_module('TOPGeometryParInitializer')
332 
333 # Time Recalibrator
334 main.add_module('TOPTimeRecalibrator', subtractBunchTime=False)
335 
336 # Channel masking
337 main.add_module('TOPChannelMasker')
338 
339 # Exclude BS13d
340 main.add_module(Mask_BS13d())
341 
342 # Bunch finder
343 main.add_module('TOPBunchFinder', usePIDLikelihoods=True, subtractRunningOffset=False)
344 
345 # Module T0 calibrator
346 main.add_module(ModuleT0cal())
347 
348 # Print progress
349 main.add_module('Progress')
350 
351 # Process events
352 b2.process(main)
353 
354 # Print statistics
355 print(b2.statistics)
Class to access a DBObjPtr from Python.
Definition: PyDBObj.h:48
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:56
A
matrix denoting which slot pair was involved
h_slotPairs_acc
histogram to store Chi2/ndf of successfully fitted time differences
h_relModuleT0
histogram to store minimization results (relative module T0)
h_moduleT0
histogram to store final results (module T0)
sigma
uncertainties of fitted time differences
h_out
list of histograms selected for output
h_deltaT0
histograms of time difference between two slots
h_slotPairs
histogram of number of entries per slot pairs