Belle II Software  release-06-02-00
beamspot.py
1 # -*- coding: utf-8 -*-
2 
3 
10 
11 '''
12 Validation of the Beam Spot calibration
13 '''
14 
15 
16 import basf2
17 from prompt import ValidationSettings
18 import ROOT
19 import sys
20 import subprocess
21 import json
22 
23 import numpy as np
24 import scipy.linalg as la
25 import scipy.stats
26 import matplotlib.pyplot as plt
27 
28 import re
29 import os
30 from glob import glob
31 from math import sqrt, frexp, asin
32 from itertools import groupby
33 
34 from datetime import datetime, timedelta
35 from ROOT.Belle2 import Unit
36 
37 
38 settings = ValidationSettings(name='BeamSpot Calibrations',
39  description=__doc__,
40  download_files=['stdout'],
41  expert_config={})
42 
43 
44 # Get the eigen parameters of the matrix
45 def getEigenPars(covM):
46  sxx, syy, szz, sxy, sxz, syz = covM
47 
48  mat = np.array([[sxx, sxy, sxz],
49  [sxy, syy, syz],
50  [sxz, syz, szz]])
51 
52  result = la.eig(mat)
53  eigVals = result[0].real
54  eigVecs = result[1]
55 
56  # sort it by size
57  idx = eigVals.argsort()[::-1]
58  eigVals = eigVals[idx]
59  eigVecs = eigVecs[:, idx]
60 
61  assert(eigVals[0] > 1e-4)
62  assert(eigVals[1] > 1e-4)
63  assert(eigVals[2] > 1e-4)
64 
65  sxEig = sqrt((eigVals[1]))
66  syEig = sqrt((eigVals[2]))
67  szEig = sqrt((eigVals[0]))
68 
69  # get the signs right
70  if eigVecs[2][0] < 0:
71  eigVecs[:, 0] *= -1
72 
73  if eigVecs[0][1] < 0:
74  eigVecs[:, 1] *= -1
75 
76  if eigVecs[1][2] < 0:
77  eigVecs[:, 2] *= -1
78 
79  # calculate the angles in mrad
80  angleXZ = Unit.rad / Unit.mrad * asin(eigVecs[0][0] / sqrt(1 - eigVecs[1][0]**2))
81  angleYZ = Unit.rad / Unit.mrad * asin(eigVecs[1][0])
82  angleXY = Unit.rad / Unit.mrad * asin(eigVecs[1][1] / sqrt(1 - eigVecs[1][0]**2))
83 
84  return [sxEig, syEig, szEig, angleXZ, angleYZ, angleXY]
85 
86 
87 # Load the values of the Beam Spot parameters into the list
88 def getBSvalues(path):
89  runDict = {}
90  with open(path + '/database.txt') as fDB:
91  for ll in fDB:
92  ll = ll.strip()
93  ll = ll.split(' ')
94 
95  # deriving file name
96  fN = ll[0].replace('/', '_') + '_rev_' + ll[1] + '.root'
97 
98  r = ll[2].split(',')
99  runDict[fN] = (int(r[0]), int(r[1]))
100 
101  arr = []
102 
103  fList = glob(path + '/*.root')
104  for fName in fList:
105  bName = os.path.basename(fName)
106  runExp = runDict[bName]
107 
108  f = ROOT.TFile.Open(fName)
109  bsAll = f.Get("BeamSpot")
110  assert(bsAll.ClassName() == "Belle2::EventDependency")
111 
112  evNums = bsAll.getEventNumbers()
113 
114  for i in range(len(evNums) + 1):
115  bs = bsAll.getObjectByIndex(i)
116  ipR = bs.getIPPosition()
117  ip = [ipR(i) * Unit.cm / Unit.um for i in range(3)] # from cm to um
118  ipeR = bs.getIPPositionCovMatrix()
119  ipe = [sqrt(ipeR(i, i)) * Unit.cm / Unit.um for i in range(3)] # from cm to um
120  covR = bs.getSizeCovMatrix()
121  covM = (covR(0, 0), covR(1, 1), covR(2, 2), covR(0, 1), covR(0, 2), covR(1, 2))
122  covM = [x * (Unit.cm / Unit.um)**2 for x in covM] # from cm2 to um2
123 
124  tStart = ipeR(0, 1) * 1e20
125  tEnd = ipeR(0, 2) * 1e20
126  tStart = int(frexp(ipeR(0, 1))[0] * 2**53) % 2**32 / 3600.
127  tEnd = int(frexp(ipeR(0, 2))[0] * 2**53) % 2**32 / 3600.
128 
129  eigM = getEigenPars(covM)
130  arr.append((runExp, tStart, tEnd, ip, ipe, covM, eigM))
131  f.Close()
132 
133  arr = sorted(arr, key=lambda x: x[1])
134 
135  return arr
136 
137 # Print the Beam Spot parameters to the text file
138 
139 
140 def printToFile(arr):
141 
142  with open('bsData.csv', 'w') as outFile:
143  outFile.write('exp run tStart tEnd xIP yIP zIP xeIP yeIP zeIP sXX sYY sZZ sXY sXZ '
144  + 'sYZ xSizeEig ySizeEig zSizeEig xzAngle yzAngle xyAngle\n')
145 
146  for e in arr:
147  outFile.write(
148  # exp run
149  str(e[0][0]) + ' ' +
150  str(e[0][1]) + ' ' +
151  # tStart tEnd
152  str(e[1]) + ' ' +
153  str(e[2]) + ' ' +
154  # xIP yIP zIP
155  str(e[3][0]) + ' ' +
156  str(e[3][1]) + ' ' +
157  str(e[3][2]) + ' ' +
158  # xeIP yeIP zeIP
159  str(e[4][0]) + ' ' +
160  str(e[4][1]) + ' ' +
161  str(e[4][2]) + ' ' +
162  # sXX sYY sZZ sXY sXZ sYZ
163  str(e[5][0]) + ' ' +
164  str(e[5][1]) + ' ' +
165  str(e[5][2]) + ' ' +
166  str(e[5][3]) + ' ' +
167  str(e[5][4]) + ' ' +
168  str(e[5][5]) + ' ' +
169  # xSizeEig ySizeEig zSizeEig xzAngle yzAngle xyAngle
170  str(e[6][0]) + ' ' +
171  str(e[6][1]) + ' ' +
172  str(e[6][2]) + ' ' +
173  str(e[6][3]) + ' ' +
174  str(e[6][4]) + ' ' +
175  str(e[6][5]) + '\n')
176  return arr
177 
178 
179 # Create a plot with the specified variable
180 def plotVar(arr, limits, vName, getterV, getterE=None):
181 
182  tVals = []
183  Vals = []
184  Errs = []
185 
186  tGapVals = []
187  GapVals = []
188 
189  for i, el in enumerate(arr):
190  s, e = el[1], el[2]
191  s = datetime.utcfromtimestamp((s + 9) * 3600) # Convert to the JST (+9 hours)
192  e = datetime.utcfromtimestamp((e + 9) * 3600)
193  tVals.append(s)
194  tVals.append(e)
195  Vals.append(getterV(el))
196  Vals.append(getterV(el))
197 
198  if getterE is not None:
199  Errs.append(getterE(el))
200  Errs.append(getterE(el))
201 
202  # Add breaks for longer gap if not the last interval
203  if i >= len(arr) - 1:
204  continue
205 
206  dt = (arr[i + 1][1] - arr[i][2]) * 3600
207 
208  # only consider gaps longer than 10 mins
209  if dt < 10 * 60:
210  continue
211 
212  # start-time of gap, end-time of gap
213  gS = datetime.utcfromtimestamp((arr[i][2] + 9) * 3600) # Convert to the JST (+9 hours)
214  gE = datetime.utcfromtimestamp((arr[i + 1][1] + 9) * 3600)
215 
216  tVals.append(gS + (gE - gS) / 2)
217  Vals.append(np.nan)
218  if getterE is not None:
219  Errs.append(np.nan)
220 
221  # store curve connecting gaps
222  tGapVals.append(gS - timedelta(seconds=1))
223  tGapVals.append(gS)
224  tGapVals.append(gE)
225  tGapVals.append(gE + timedelta(seconds=1))
226 
227  GapVals.append(np.nan)
228  GapVals.append(getterV(arr[i]))
229  GapVals.append(getterV(arr[i + 1]))
230  GapVals.append(np.nan)
231 
232  plt.plot(tVals, Vals, linewidth=2, color='C0')
233  plt.plot(tGapVals, GapVals, linewidth=2, color='C0', alpha=0.35)
234 
235  Vals = np.array(Vals)
236  Errs = np.array(Errs)
237 
238  if getterE is not None:
239  plt.fill_between(tVals, Vals - Errs, Vals + Errs, alpha=0.2)
240 
241  plt.xlabel('time')
242  if 'Angle' in vName:
243  plt.ylabel(vName + ' [mrad]')
244  else:
245  plt.ylabel(vName + ' [um]')
246 
247  loc = 'plots/allData'
248  if limits is not None:
249  plt.xlim(datetime.strptime(limits[0], '%Y-%m-%d'), datetime.strptime(limits[1], '%Y-%m-%d'))
250  loc = 'plots/' + limits[0] + 'to' + limits[1]
251 
252  # if not os.path.isdir(loc):
253  os.makedirs(loc, exist_ok=True)
254 
255  plt.savefig(loc + '/' + vName + '.png')
256  plt.clf()
257 
258 
259 def plotPullSpectrum(arr, vName, getterV, getterE):
260 
261  valsOrg = np.array([getterV(v) for v in arr])
262  errsOrg = np.array([getterE(v) for v in arr])
263 
264  # remove repeating values
265  indx = np.array([list(g)[0][0] for k, g in groupby(zip(np.arange(len(valsOrg)), valsOrg), lambda x: x[1])])
266  vals = valsOrg[indx]
267  errs = errsOrg[indx]
268 
269  diffs = (vals[1:] - vals[:-1]) / np.hypot(errs[1:], errs[:-1])
270 
271  # Get 1sigma quantiles of the normal distribution
272  n1, n2 = scipy.stats.norm.cdf([-1, 1])
273 
274  q1 = np.quantile(diffs, n1)
275  q2 = np.quantile(diffs, n2)
276 
277  plt.hist(diffs, bins=np.linspace(-20, 20, 80), label=f'q1={round(q1,1)}, q2={round(q2,1)}')
278  plt.xlabel(vName)
279  plt.legend()
280 
281  loc = 'plots/allData'
282  plt.savefig(loc + '/' + vName + '.png')
283  plt.clf()
284 
285 
286 def run_validation(job_path, input_data_path, requested_iov, expert_config):
287  '''
288  Run the validation.
289  '''
290 
291  # Expert config can contain the time ranges of the plots
292  if expert_config != '':
293  expert_config = json.loads(expert_config)
294 
295  allLimits = [None]
296  if expert_config is not None and 'plotsRanges' in expert_config:
297  allLimits += expert_config['plotsRanges']
298 
299  # Path to the database.txt file and to the payloads.
300  dbFile = glob(f'{job_path}/**/database.txt', recursive=True)
301  dbFile = [db for db in dbFile if 'algorithm_output' not in db]
302  assert(len(dbFile) == 1)
303  dbFile = dbFile[0]
304  inputDir = dbFile[:dbFile.rfind('/')]
305 
306  arr = getBSvalues(inputDir)
307 
308  # print the results to the CSV file
309  printToFile(arr)
310 
311  plt.figure(figsize=(18, 9))
312 
313  # plot the results
314  for limits in allLimits:
315  plotVar(arr, limits, 'xIP', lambda e: e[3][0], lambda e: e[4][0])
316  plotVar(arr, limits, 'yIP', lambda e: e[3][1], lambda e: e[4][1])
317  plotVar(arr, limits, 'zIP', lambda e: e[3][2], lambda e: e[4][2])
318 
319  plotVar(arr, limits, 'xSpotSize', lambda e: sqrt(e[5][0]))
320  plotVar(arr, limits, 'ySpotSize', lambda e: sqrt(e[5][1]))
321  plotVar(arr, limits, 'zSpotSize', lambda e: sqrt(e[5][2]))
322 
323  plotVar(arr, limits, 'xSpotSizeEig', lambda e: e[6][0])
324  plotVar(arr, limits, 'ySpotSizeEig', lambda e: e[6][1])
325  plotVar(arr, limits, 'zSpotSizeEig', lambda e: e[6][2])
326 
327  plotVar(arr, limits, 'xzSpotAngle', lambda e: e[6][3])
328  plotVar(arr, limits, 'yzSpotAngle', lambda e: e[6][4])
329  plotVar(arr, limits, 'xySpotAngle', lambda e: e[6][5])
330 
331  plotPullSpectrum(arr, 'xIPpulls', lambda e: e[3][0], lambda e: e[4][0])
332  plotPullSpectrum(arr, 'yIPpulls', lambda e: e[3][1], lambda e: e[4][1])
333  plotPullSpectrum(arr, 'zIPpulls', lambda e: e[3][2], lambda e: e[4][2])
334 
335 
336 if __name__ == "__main__":
337  run_validation(*sys.argv[1:])