Belle II Software development
beamparameters.py
1#!/usr/bin/env python3
2
3
10
11"""
12Scripts and functions to set the BeamParameters from known configurations
13"""
14
15from basf2 import B2FATAL, B2INFO
16import math
17import numpy as np
18
19
20beamparameter_defaults = {
21 "energyHER": 5.28,
22 "energyLER": 5.28,
23 "angleXHER": 0,
24 "angleXLER": math.pi,
25 "covHER": [0],
26 "covLER": [0],
27 "vertex": [0, 0, 0],
28 "covVertex": [0]
29}
30
31
39beamparameter_presets = {
40 "SuperKEKB": (None, {
41 # beam energies in GeV
42 "energyHER": 7.,
43 "energyLER": 4.,
44 # beam angles with respect to the z axis in radian, negative gets
45 # converted into pi - |angle|
46 "angleXHER": 0.0415,
47 "angleXLER": -0.0415,
48 # covariance matrices for the beam parametrized in
49 # (E, theta_x, theta_y) = energy and horizontal and vertical angular
50 # spread
51 "covHER": np.array([0.00513]) ** 2,
52 "covLER": np.array([0.002375]) ** 2,
53 # bunch size in cm
54 "bunchHER": np.array([10.2e-3, 59e-6, 6.0]) * 1e-1,
55 "bunchLER": np.array([7.75e-3, 59e-6, 5.0]) * 1e-1,
56 }),
57 "Y1S": ("SuperKEKB", { # m(Y1S) = 9.460 GeV
58 "energyHER": 6.263,
59 "energyLER": 3.579,
60 }),
61 "Y2S": ("SuperKEKB", { # m(Y2S) = 10.023 GeV
62 "energyHER": 6.635,
63 "energyLER": 3.792,
64 }),
65 "Y3S": ("SuperKEKB", { # m(Y3S) = 10.355 GeV
66 "energyHER": 6.855,
67 "energyLER": 3.917,
68 }),
69 "Y4S": ("SuperKEKB", { # m(Y4S) = 10.579 GeV
70 "energyHER": 7.004,
71 "energyLER": 4.002,
72 }),
73 "Y5S": ("SuperKEKB", { # m(Y5S) = 10.876 GeV
74 "energyHER": 7.200,
75 "energyLER": 4.114,
76 }),
77 "Y1S-off": ("Y1S", { # m(Y1S) - 60 MeV = 9.400 GeV
78 "energyHER": 6.223,
79 "energyLER": 3.556,
80 }),
81 "Y2S-off": ("Y2S", { # m(Y2S) - 60 MeV = 9.963 GeV
82 "energyHER": 6.596,
83 "energyLER": 3.769,
84 }),
85 "Y3S-off": ("Y3S", { # m(Y3S) - 60 MeV = 10.295 GeV
86 "energyHER": 6.816,
87 "energyLER": 3.895,
88 }),
89 "Y4S-off": ("Y4S", { # m(Y4S) - 60 MeV = 10.519 GeV
90 "energyHER": 6.964,
91 "energyLER": 3.979,
92 }),
93 "Y5S-off": ("Y5S", { # m(Y5S) - 60 MeV = 10.816 GeV
94 "energyHER": 7.160,
95 "energyLER": 4.092,
96 }),
97 "Y1D1": ("SuperKEKB", { # m(Y1D) = 10.152 GeV
98 "energyHER": 6.720,
99 "energyLER": 3.840,
100 }),
101 "Y2D1": ("SuperKEKB", { # m(Y2D) = 10.425 GeV
102 "energyHER": 6.901,
103 "energyLER": 3.944,
104 }),
105 "Y6S": ("SuperKEKB", { # m(Y6S) = 11.020 GeV
106 "energyHER": 7.295,
107 "energyLER": 4.169,
108 }),
109 "Yb10775": ("SuperKEKB", { # m(Yb10775) = 10.775 GeV
110 "energyHER": 7.133,
111 "energyLER": 4.076,
112 }),
113 "KEKB-Belle": (None, {
114 "energyHER": 7.998213,
115 "energyLER": 3.499218,
116 "angleXHER": 0.022,
117 }),
118 "LEP1": (None, {
119 "energyHER": 45.6,
120 "energyLER": 45.6,
121 }),
122 "DAPHNE": (None, {
123 "energyHER": 1.02 / 2.,
124 "energyLER": 1.02 / 2.,
125 }),
126}
127
128
129def __get_4vector(energy, angle):
130 """Calculate an 4vector for electron/positron from energy and angle"""
131 import ROOT
132 m = 0.511e-3
133 pz = (energy ** 2 - m ** 2) ** .5
134 v = ROOT.Math.PxPyPzEVector(0, 0, pz, energy)
135 if angle < 0:
136 angle = math.pi + angle
137 rotationY = ROOT.Math.RotationY(angle)
138 v = rotationY(v)
139 return v
140
141
142def rot_matrix_y(angle):
143 """Return rotation matrix for a rotation of angle around the y-axis"""
144 c = math.cos(angle)
145 s = math.sin(angle)
146 return np.matrix([(c, 0, s), (0, 1, 0), (-s, 0, c)])
147
148
149def cov_matrix(beamsize, angle):
150 """Return the covariance matrix for the given bunch size and angle with
151 respect to the z axis"""
152 if angle < 0:
153 angle = math.pi + angle
154 rot = rot_matrix_y(angle)
155 cov = np.matrix(np.diag(beamsize ** 2))
156 return rot * cov * rot.T
157
158
159def calculate_beamspot(pos_her, pos_ler, size_her, size_ler, angle_her, angle_ler):
160 """
161 Function to calculate position and covariance matrix of the beamspot.
162
163 This function calculates the mean position and the covariance matrix for the
164 beamspot from the bunch positions, bunch sizes and the angles of the beam
165 with respect to the z-axis. It assumes normal distributed bunch densities and
166 multiplies the PDFs for the two bunches to get a resulting beamspot
167
168 negative angles are interpreted as "math.pi - abs(angle)"
169
170 Args:
171 pos_her (list): mean position of the HER bunch (in cm)
172 pos_ler (list): mean position of the LER bunch (in cm)
173 size_her (list): HER bunch size in cm
174 size_ler (list): LER bunch size in cm
175 angle_her (float): angle between HER direction and z-axis
176 angle_her (float): angle between LER direction and z-axis
177
178 Returns:
179 pos (array): array (x, y, z) with the spot center in cm
180 cov (matrix): covariance matrix of the beam spot
181 """
182
183 cov_her = cov_matrix(size_her, angle_her)
184 cov_ler = cov_matrix(size_ler, angle_ler)
185 # calculate the inverse of the sum of the two covariance matrices
186 cov_sum_i = (cov_her + cov_ler).I
187 # and use it to calculate the product of the two gaussian distributions,
188 # covariance and mean
189 cov_spot = cov_her * cov_sum_i * cov_ler
190 beampos_spot = np.array(
191 cov_ler * cov_sum_i * np.matrix(pos_her).T +
192 cov_her * cov_sum_i * np.matrix(pos_ler).T
193 ).T[0]
194 return beampos_spot, cov_spot
195
196
197def add_beamparameters(path, name, E_cms=None, **argk):
198 """Add BeamParameter module to a given path.
199 Note that if you append this to the path, you will not
200 read the configuration of BeamParameters from the database
201 Beware that this function does not mimic at all the behaviour
202 of the official MC samples.
203 Args:
204 path (basf2.Path instance): path to add the module to
205 name (str): name of the beamparameter settings to use
206 E_cms (float): center of mass energy. If not None the beamenergies will
207 be scaled accordingly to achieve E_cms
208
209 Additional keyword arguments will be passed directly to the module as parameters.
210 """
211 # calculate all the parameter values from the preset
212 values = calculate_beamparameters(name, E_cms)
213 # add a BeamParametes module to the path
214 module = path.add_module("BeamParameters")
215 module.set_name(f"BeamParameters:{name}")
216 # finally, set all parameters and return the module
217 module.param(values)
218 # and override parameters with any additional keyword arguments
219 module.param(argk)
220 return module
221
222
223def calculate_beamparameters(name, E_cms=None):
224 """Get/calculate the necessary beamparameter values from the given preset name,
225 optionally scale it to achieve the given cms energy
226
227 Args:
228 name (str): name of the beamparameter settings to use
229 E_cms (float): center of mass energy. If not None the beamenergies will
230 be scaled accordingly to achieve E_cms
231 """
232 if name not in beamparameter_presets:
233 B2FATAL("Unknown beamparameter preset: '%s', use one of %s" %
234 (name, ", ".join(sorted(beamparameter_presets.keys()))))
235 return None
236
237 # copy the defaults
238 values = beamparameter_defaults.copy()
239
240 # collect all the values from the preset by also looping over the
241 preset = {}
242
243 def collect_values(preset, name):
244 parent, values = beamparameter_presets[name]
245 if parent is not None:
246 collect_values(preset, parent)
247 preset.update(values)
248
249 collect_values(preset, name)
250
251 # and update with the preset
252 for k in values.keys():
253 if k in preset:
254 if isinstance(preset[k], np.ndarray):
255 values[k] = list(preset[k].flat)
256 else:
257 values[k] = preset[k]
258
259 # if we have bunch sizes we want to compute the vertex covariance from them
260 if "bunchHER" in preset and "bunchLER" in preset:
261 her_size = preset["bunchHER"]
262 ler_size = preset["bunchLER"]
263 her_angle = values["angleXHER"]
264 ler_angle = values["angleXLER"]
265 _, cov = calculate_beamspot(
266 [0, 0, 0], [0, 0, 0], her_size, ler_size, her_angle, ler_angle
267 )
268 values["covVertex"] = list(cov.flat)
269
270 if E_cms is not None:
271 ler = __get_4vector(values["energyLER"], values["angleXLER"])
272 her = __get_4vector(values["energyHER"], values["angleXHER"])
273 mass = (her + ler).M()
274 scale = E_cms / mass
275 B2INFO(f"Scaling beam energies by {scale:g} to obtain E_cms = {E_cms:g} GeV")
276 values["energyHER"] *= scale
277 values["energyLER"] *= scale
278
279 return values
280
281
282def _get_collisions_invariant_mass(experiment, run, verbose=False):
283 """
284 Convenient function for getting the collisions invariant mass by querying the BeamParameters object in the Conditions Database.
285 Calling this function directly might have undesired side effects: users should always use
286 :func:`get_collisions_invariant_mass()` in their steering files.
287
288 Args:
289 experiment (int): experiment number of the basf2 process
290 run (int): run number of the basf2 process
291 verbose (bool): if True, it prints some INFO messages, e.g. the invariant mass retrieved from the CDB
292 """
293 from ROOT import Belle2 as B2 # noqa
294 # Initialize the datastore and create an EventMetaData object
295 B2.DataStore.Instance().setInitializeActive(True)
296 event_meta_data = B2.PyStoreObj('EventMetaData')
297 event_meta_data.registerInDataStore()
298 event_meta_data.assign(B2.EventMetaData(0, run, experiment), True)
299 B2.DataStore.Instance().setInitializeActive(False)
300 # Get the invariant mass by querying the database
301 beams = B2.PyDBObj('BeamParameters')
302 if not beams.isValid():
303 B2FATAL('BeamParameters is not valid')
304 collisions_invariant_mass = beams.getMass()
305 # Print some messages if requested
306 if verbose:
307 B2INFO('Info for BeamParameters:\n'
308 f' collisions invariant mass [GeV]: {collisions_invariant_mass}\n'
309 f' IoV: {beams.getIoV().getExperimentLow()}, {beams.getIoV().getRunLow()}, '
310 f'{beams.getIoV().getExperimentHigh()}. {beams.getIoV().getRunHigh()}\n'
311 f' revision: {beams.getRevision()}\n'
312 f' globaltag: {beams.getGlobaltag()}')
313 # Before returning, run some cleanup
314 B2.DataStore.Instance().reset()
315 B2.Database.reset()
316 return collisions_invariant_mass
317
318
319def get_collisions_invariant_mass(experiment, run, verbose=False):
320 """
321 Convenient function for getting the collisions invariant mass by querying the BeamParameters object in the Conditions Database.
322 It is useful for setting the center of mass energy when producing Monte Carlo samples with external generators (e.g. MadGraph,
323 WHIZARD) which are not directly interfaced with basf2.
324
325 Args:
326 experiment (int): experiment number of the basf2 process
327 run (int): run number of the basf2 process
328 verbose (bool): if True, it prints some INFO messages, e.g. the invariant mass retrieved from the CDB
329 """
330
331 def _run_function_in_process(func, *args, **kwargs):
332 """
333 Small helper function for running another function in a process and thus avoiding undesired side effects.
334 """
335 from multiprocessing import Pool # noqa
336 # Create a Pool with a single process
337 with Pool(processes=1) as pool:
338 # Use apply_async to execute the function and get the result
339 result = pool.apply_async(func, args, kwargs)
340 # Wait for the result and retrieve it
341 output = result.get()
342 return output
343
344 # Run _get_collisions_invariant_mass with the helper function to avoid issues with basf2
345 return _run_function_in_process(
346 func=_get_collisions_invariant_mass, experiment=experiment, run=run, verbose=verbose
347 )
348
349
350if __name__ == "__main__":
351 # if called directly we will
352 # 1. calculate the beamspot for the SuperKEKB preset and display it if possible
353 # 2. calculate beam energies for Y1S - Y5S and offresonance
354 np.set_printoptions(precision=3)
355
356 # calculate beamspot for SuperKEKB for testing
357
358 values = beamparameter_presets["SuperKEKB"][1]
359 beampos_spot, cov_spot = calculate_beamspot(
360 [0, 0, 0], [0, 0, 0],
361 values["bunchHER"], values["bunchLER"],
362 values["angleXHER"], values["angleXLER"],
363 )
364 print("Beamspot position:")
365 print(beampos_spot)
366 print("Beamspot covariance:")
367 print(cov_spot)
368 print("Beamspot dimensions (in mm, axes in arbitrary order):")
369 print(np.linalg.eig(cov_spot)[0] ** .5 * 10)
370
371 # see if we can plot it
372 try:
373 from matplotlib import pyplot as pl
374
375 cov_her = cov_matrix(values["bunchHER"], values["angleXHER"])
376
377 cov_ler = cov_matrix(values["bunchLER"], values["angleXLER"])
378
379 points_her = np.random.multivariate_normal([0, 0, 0], cov_her, 1000)
380
381 points_ler = np.random.multivariate_normal([0, 0, 0], cov_ler, 1000)
382
383 points_spot = np.random.multivariate_normal(beampos_spot, cov_spot, 1000)
384 pl.scatter(points_her[:, 2], points_her[:, 0], s=5, marker=".", c="b", edgecolors="None", label="HER")
385 pl.scatter(points_ler[:, 2], points_ler[:, 0], s=5, marker=".", c="r", edgecolors="None", label="LER")
386 pl.scatter(points_spot[:, 2], points_spot[:, 0], s=5, marker=".", c="g", edgecolors="None", label="spot")
387 pl.legend()
388 pl.xlabel("z / cm")
389 pl.ylabel("x / cm")
390 pl.tight_layout()
391 pl.show()
392 except ImportError:
393 pass
394
395
396 eher = values["energyHER"]
397
398 eler = values["energyLER"]
399
400 aher = values["angleXHER"]
401
402 aler = values["angleXLER"]
403
404 # calculate beam energies for Y1S - Y4S by scaling them from the nominal
405 # beam energies for the SuperKEKB preset
406
407 ler = __get_4vector(eher, aher)
408
409 her = __get_4vector(eler, aler)
410
411 mass = (her + ler).M()
412 print("Nominal CMS Energy: ", mass)
413
414
415 targets = {
416 "Y1S": 9460.30e-3,
417 "Y2S": 10023.26e-3,
418 "Y3S": 10355.2e-3,
419 "Y4S": 10579.4e-3,
420 "Y5S": 10876e-3
421 }
422 for name, energy in sorted(targets.items()):
423
424 scale = energy / mass
425 print(f""" "{name}": ("SuperKEKB", {{ # m({name}) = {energy:.3f} GeV
426 "energyHER": {eher * scale:.3f},
427 "energyLER": {eler * scale:.3f},
428 }}),""")
429
430 for name, energy in sorted(targets.items()):
431
432 scale = (energy - 60e-3) / mass
433 print("""\
434 "{}-off": ("{}", {{ # m({}) - 60 MeV = {:.3f} GeV
435 "energyHER": {:.3f},
436 "energyLER": {:.3f},
437 }}),""".format(name, name, name, energy - 60e-3, eher * scale, eler * scale))