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