Belle II Software  release-06-02-00
boostvector.py
1 # -*- coding: utf-8 -*-
2 
3 
10 
11 '''
12 Validation of the Boost Vector 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, hypot, atan2
32 
33 from datetime import datetime, timedelta
34 from ROOT.Belle2 import Unit
35 
36 
37 settings = ValidationSettings(name='BoostVector Calibrations',
38  description=__doc__,
39  download_files=['stdout'],
40  expert_config={})
41 
42 
43 # Get the eigen parameters of the matrix
44 def getBVvalsErrs(bv, bve):
45 
46  bv = np.array(bv)
47  bve = np.array(bve)
48 
49  mag = np.linalg.norm(bv)
50  xzAngle = atan2(bv[0], bv[2])
51  yzAngle = atan2(bv[1], bv[2])
52 
53  magUnc = np.linalg.norm(bv * bve) / mag
54  xzAngleUnc = 1. / (1 + (bv[0] / bv[2])**2) * hypot(bve[0] / bv[2], -bve[1] * bv[0] / bv[2]**2)
55  yzAngleUnc = 1. / (1 + (bv[1] / bv[2])**2) * hypot(bve[1] / bv[2], -bve[1] * bv[1] / bv[2]**2)
56 
57  tomrad = Unit.rad / Unit.mrad
58  return [1e3 * mag, tomrad * xzAngle, tomrad * yzAngle], [1e3 * magUnc, tomrad * xzAngleUnc, tomrad * yzAngleUnc]
59 
60 
61 # Load the values of the Boost Vector parameters into the list
62 def getBVvalues(path):
63  runDict = {}
64  with open(path + '/database.txt') as fDB:
65  for ll in fDB:
66  ll = ll.strip()
67  ll = ll.split(' ')
68 
69  # deriving file name
70  fN = ll[0].replace('/', '_') + '_rev_' + ll[1] + '.root'
71 
72  r = ll[2].split(',')
73  runDict[fN] = (int(r[0]), int(r[1]))
74 
75  arr = []
76 
77  fList = glob(path + '/*.root')
78  for fName in fList:
79  bName = os.path.basename(fName)
80  runExp = runDict[bName]
81 
82  f = ROOT.TFile.Open(fName)
83  bvAll = f.Get("CollisionBoostVector")
84  assert(bvAll.ClassName() == "Belle2::EventDependency")
85 
86  evNums = bvAll.getEventNumbers()
87 
88  for i in range(len(evNums) + 1):
89  bvObj = bvAll.getObjectByIndex(i)
90 
91  bvV = bvObj.getBoost()
92  bv = [bvV(i) for i in range(3)]
93  bveV = bvObj.getBoostCovariance()
94  bve = [sqrt(bveV(i, i)) for i in range(3)]
95 
96  tStart = int(frexp(bveV(0, 1))[0] * 2**53) % 2**32 / 3600.
97  tEnd = int(frexp(bveV(0, 2))[0] * 2**53) % 2**32 / 3600.
98 
99  bvV, bvE = getBVvalsErrs(bv, bve)
100 
101  arr.append((runExp, tStart, tEnd, bvV, bvE))
102  f.Close()
103 
104  arr = sorted(arr, key=lambda x: x[1])
105 
106  return arr
107 
108 # Print the Boost Vector parameters to the text file
109 
110 
111 def printToFile(arr):
112 
113  with open('bvData.csv', 'w') as outFile:
114  outFile.write('exp run tStart tEnd mag xzAngle yzAngle magUnc xzAngleUnc yzAngleUnc \n')
115 
116  for e in arr:
117  outFile.write(
118  # exp run
119  str(e[0][0]) + ' ' +
120  str(e[0][1]) + ' ' +
121  # tStart tEnd
122  str(e[1]) + ' ' +
123  str(e[2]) + ' ' +
124  # mag xzAngle yzAngle
125  str(e[3][0]) + ' ' +
126  str(e[3][1]) + ' ' +
127  str(e[3][2]) + ' ' +
128  # magUnc xzAngleUnc yzAngleUnc
129  str(e[4][0]) + ' ' +
130  str(e[4][1]) + ' ' +
131  str(e[4][2]) + '\n')
132  return arr
133 
134 
135 # Create a plot with the specified variable
136 def plotVar(arr, limits, vName, getterV, getterE=None):
137 
138  tVals = []
139  Vals = []
140  Errs = []
141 
142  tGapVals = []
143  GapVals = []
144 
145  for i, el in enumerate(arr):
146  s, e = el[1], el[2]
147  s = datetime.utcfromtimestamp((s + 9) * 3600) # Convert to the JST (+9 hours)
148  e = datetime.utcfromtimestamp((e + 9) * 3600)
149  tVals.append(s)
150  tVals.append(e)
151  Vals.append(getterV(el))
152  Vals.append(getterV(el))
153 
154  if getterE is not None:
155  Errs.append(getterE(el))
156  Errs.append(getterE(el))
157 
158  # Add breaks for longer gap if not the last interval
159  if i >= len(arr) - 1:
160  continue
161 
162  dt = (arr[i + 1][1] - arr[i][2]) * 3600
163 
164  # only consider gaps longer than 10 mins
165  if dt < 10 * 60:
166  continue
167 
168  # start-time of gap, end-time of gap
169  gS = datetime.utcfromtimestamp((arr[i][2] + 9) * 3600) # Convert to the JST (+9 hours)
170  gE = datetime.utcfromtimestamp((arr[i + 1][1] + 9) * 3600)
171 
172  tVals.append(gS + (gE - gS) / 2)
173  Vals.append(np.nan)
174  if getterE is not None:
175  Errs.append(np.nan)
176 
177  # store curve connecting gaps
178  tGapVals.append(gS - timedelta(seconds=1))
179  tGapVals.append(gS)
180  tGapVals.append(gE)
181  tGapVals.append(gE + timedelta(seconds=1))
182 
183  GapVals.append(np.nan)
184  GapVals.append(getterV(arr[i]))
185  GapVals.append(getterV(arr[i + 1]))
186  GapVals.append(np.nan)
187 
188  plt.plot(tVals, Vals, linewidth=2, color='C0')
189  plt.plot(tGapVals, GapVals, linewidth=2, color='C0', alpha=0.35)
190 
191  Vals = np.array(Vals)
192  Errs = np.array(Errs)
193 
194  if getterE is not None:
195  plt.fill_between(tVals, Vals - Errs, Vals + Errs, alpha=0.2)
196 
197  plt.xlabel('time')
198  if 'Angle' in vName:
199  plt.ylabel(vName + ' [mrad]')
200  else:
201  plt.ylabel(vName + r' $[10^{-3}]$')
202 
203  loc = 'plots/allData'
204  if limits is not None:
205  plt.xlim(datetime.strptime(limits[0], '%Y-%m-%d'), datetime.strptime(limits[1], '%Y-%m-%d'))
206  loc = 'plots/' + limits[0] + 'to' + limits[1]
207 
208  # if not os.path.isdir(loc):
209  os.makedirs(loc, exist_ok=True)
210 
211  plt.savefig(loc + '/' + vName + '.png')
212  plt.clf()
213 
214 
215 def plotPullSpectrum(arr, vName, getterV, getterE):
216 
217  from itertools import groupby
218  vals = np.array([k for k, g in groupby([getterV(v) for v in arr])])
219  errs = np.array([k for k, g in groupby([getterE(v) for v in arr])])
220  assert(len(vals) == len(errs))
221 
222  diffs = (vals[1:] - vals[:-1]) / np.hypot(errs[1:], errs[:-1])
223 
224  # Get 1sigma quantiles of the normal distribution
225  n1, n2 = scipy.stats.norm.cdf([-1, 1])
226 
227  q1 = np.quantile(diffs, n1)
228  q2 = np.quantile(diffs, n2)
229 
230  plt.hist(diffs, bins=np.linspace(-10, 10, 20), label=f'q1={round(q1,1)}, q2={round(q2,1)}')
231  plt.xlabel(vName)
232  plt.legend()
233 
234  loc = 'plots/allData'
235  plt.savefig(loc + '/' + vName + '.png')
236  plt.clf()
237 
238 
239 def run_validation(job_path, input_data_path, requested_iov, expert_config):
240  '''
241  Run the validation.
242  '''
243 
244  # Expert config can contain the time ranges of the plots
245  if expert_config != '':
246  expert_config = json.loads(expert_config)
247 
248  allLimits = [None]
249  if expert_config is not None and 'plotsRanges' in expert_config:
250  allLimits += expert_config['plotsRanges']
251 
252  # Path to the database.txt file and to the payloads.
253  dbFile = glob(f'{job_path}/**/database.txt', recursive=True)
254  dbFile = [db for db in dbFile if 'algorithm_output' not in db]
255  assert(len(dbFile) == 1)
256  dbFile = dbFile[0]
257  inputDir = dbFile[:dbFile.rfind('/')]
258 
259  arr = getBVvalues(inputDir)
260 
261  # print the results to the CSV file
262  printToFile(arr)
263 
264  plt.figure(figsize=(18, 9))
265 
266  # plot the results
267  for limits in allLimits:
268  plotVar(arr, limits, 'boostMag', lambda e: e[3][0], lambda e: e[4][0])
269  plotVar(arr, limits, 'xzBoostAngle', lambda e: e[3][1], lambda e: e[4][1])
270  plotVar(arr, limits, 'yzBoostAngle', lambda e: e[3][2], lambda e: e[4][2])
271 
272  plotPullSpectrum(arr, 'boostMagPulls', lambda e: e[3][0], lambda e: e[4][0])
273  plotPullSpectrum(arr, 'xzBoostAnglePulls', lambda e: e[3][1], lambda e: e[4][1])
274  plotPullSpectrum(arr, 'yzBoostAnglePulls', lambda e: e[3][2], lambda e: e[4][2])
275 
276 
277 if __name__ == "__main__":
278  run_validation(*sys.argv[1:])