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