12 Validation of the Beam Spot calibration
17 from prompt
import ValidationSettings
24 import scipy.linalg
as la
26 import matplotlib.pyplot
as plt
31 from math
import sqrt, frexp, asin
32 from itertools
import groupby
34 from datetime
import datetime, timedelta
35 from ROOT.Belle2
import Unit
38 settings = ValidationSettings(name=
'BeamSpot Calibrations',
40 download_files=[
'stdout'],
45 def getEigenPars(covM):
46 sxx, syy, szz, sxy, sxz, syz = covM
48 mat = np.array([[sxx, sxy, sxz],
53 eigVals = result[0].real
57 idx = eigVals.argsort()[::-1]
58 eigVals = eigVals[idx]
59 eigVecs = eigVecs[:, idx]
61 assert(eigVals[0] > 1e-4)
62 assert(eigVals[1] > 1e-4)
63 assert(eigVals[2] > 1e-4)
65 sxEig = sqrt((eigVals[1]))
66 syEig = sqrt((eigVals[2]))
67 szEig = sqrt((eigVals[0]))
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))
84 return [sxEig, syEig, szEig, angleXZ, angleYZ, angleXY]
88 def getBSvalues(path):
90 with open(path +
'/database.txt')
as fDB:
96 fN = ll[0].replace(
'/',
'_') +
'_rev_' + ll[1] +
'.root'
99 runDict[fN] = (int(r[0]), int(r[1]))
103 fList = glob(path +
'/*.root')
105 bName = os.path.basename(fName)
106 runExp = runDict[bName]
108 f = ROOT.TFile.Open(fName)
109 bsAll = f.Get(
"BeamSpot")
110 assert(bsAll.ClassName() ==
"Belle2::EventDependency")
112 evNums = bsAll.getEventNumbers()
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)]
118 ipeR = bs.getIPPositionCovMatrix()
119 ipe = [sqrt(ipeR(i, i)) * Unit.cm / Unit.um
for i
in range(3)]
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]
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.
129 eigM = getEigenPars(covM)
130 arr.append((runExp, tStart, tEnd, ip, ipe, covM, eigM))
133 arr = sorted(arr, key=
lambda x: x[1])
140 def printToFile(arr):
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')
180 def plotVar(arr, limits, vName, getterV, getterE=None):
189 for i, el
in enumerate(arr):
191 s = datetime.utcfromtimestamp((s + 9) * 3600)
192 e = datetime.utcfromtimestamp((e + 9) * 3600)
195 Vals.append(getterV(el))
196 Vals.append(getterV(el))
198 if getterE
is not None:
199 Errs.append(getterE(el))
200 Errs.append(getterE(el))
203 if i >= len(arr) - 1:
206 dt = (arr[i + 1][1] - arr[i][2]) * 3600
213 gS = datetime.utcfromtimestamp((arr[i][2] + 9) * 3600)
214 gE = datetime.utcfromtimestamp((arr[i + 1][1] + 9) * 3600)
216 tVals.append(gS + (gE - gS) / 2)
218 if getterE
is not None:
222 tGapVals.append(gS - timedelta(seconds=1))
225 tGapVals.append(gE + timedelta(seconds=1))
227 GapVals.append(np.nan)
228 GapVals.append(getterV(arr[i]))
229 GapVals.append(getterV(arr[i + 1]))
230 GapVals.append(np.nan)
232 plt.plot(tVals, Vals, linewidth=2, color=
'C0')
233 plt.plot(tGapVals, GapVals, linewidth=2, color=
'C0', alpha=0.35)
235 Vals = np.array(Vals)
236 Errs = np.array(Errs)
238 if getterE
is not None:
239 plt.fill_between(tVals, Vals - Errs, Vals + Errs, alpha=0.2)
243 plt.ylabel(vName +
' [mrad]')
245 plt.ylabel(vName +
' [um]')
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]
253 os.makedirs(loc, exist_ok=
True)
255 plt.savefig(loc +
'/' + vName +
'.png')
259 def plotPullSpectrum(arr, vName, getterV, getterE):
261 valsOrg = np.array([getterV(v)
for v
in arr])
262 errsOrg = np.array([getterE(v)
for v
in arr])
265 indx = np.array([list(g)[0][0]
for k, g
in groupby(zip(np.arange(len(valsOrg)), valsOrg),
lambda x: x[1])])
269 diffs = (vals[1:] - vals[:-1]) / np.hypot(errs[1:], errs[:-1])
272 n1, n2 = scipy.stats.norm.cdf([-1, 1])
274 q1 = np.quantile(diffs, n1)
275 q2 = np.quantile(diffs, n2)
277 plt.hist(diffs, bins=np.linspace(-20, 20, 80), label=f
'q1={round(q1,1)}, q2={round(q2,1)}')
281 loc =
'plots/allData'
282 plt.savefig(loc +
'/' + vName +
'.png')
286 def run_validation(job_path, input_data_path, requested_iov, expert_config):
292 if expert_config !=
'':
293 expert_config = json.loads(expert_config)
296 if expert_config
is not None and 'plotsRanges' in expert_config:
297 allLimits += expert_config[
'plotsRanges']
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)
304 inputDir = dbFile[:dbFile.rfind(
'/')]
306 arr = getBSvalues(inputDir)
311 plt.figure(figsize=(18, 9))
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])
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]))
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])
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])
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])
336 if __name__ ==
"__main__":
337 run_validation(*sys.argv[1:])