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