Belle II Software  release-08-00-10
beamparameters.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 """
5 Scripts and functions to set the BeamParameters from known configurations
6 """
7 
8 from basf2 import B2FATAL, B2INFO
9 import math
10 import numpy as np
11 
12 
13 beamparameter_defaults = {
14  "energyHER": 5.28,
15  "energyLER": 5.28,
16  "angleXHER": 0,
17  "angleXLER": math.pi,
18  "covHER": [0],
19  "covLER": [0],
20  "vertex": [0, 0, 0],
21  "covVertex": [0]
22 }
23 
24 
32 beamparameter_presets = {
33  "SuperKEKB": (None, {
34  # beam energies in GeV
35  "energyHER": 7.,
36  "energyLER": 4.,
37  # beam angles with respect to the z axis in radian, negative gets
38  # converted into pi - |angle|
39  "angleXHER": 0.0415,
40  "angleXLER": -0.0415,
41  # covariance matrices for the beam parametrized in
42  # (E, theta_x, theta_y) = energy and horizontal and vertical angular
43  # spread
44  "covHER": np.array([0.00513]) ** 2,
45  "covLER": np.array([0.002375]) ** 2,
46  # bunch size in cm
47  "bunchHER": np.array([10.2e-3, 59e-6, 6.0]) * 1e-1,
48  "bunchLER": np.array([7.75e-3, 59e-6, 5.0]) * 1e-1,
49  }),
50  "Y1S": ("SuperKEKB", { # m(Y1S) = 9.460 GeV
51  "energyHER": 6.263,
52  "energyLER": 3.579,
53  }),
54  "Y2S": ("SuperKEKB", { # m(Y2S) = 10.023 GeV
55  "energyHER": 6.635,
56  "energyLER": 3.792,
57  }),
58  "Y3S": ("SuperKEKB", { # m(Y3S) = 10.355 GeV
59  "energyHER": 6.855,
60  "energyLER": 3.917,
61  }),
62  "Y4S": ("SuperKEKB", { # m(Y4S) = 10.579 GeV
63  "energyHER": 7.004,
64  "energyLER": 4.002,
65  }),
66  "Y5S": ("SuperKEKB", { # m(Y5S) = 10.876 GeV
67  "energyHER": 7.200,
68  "energyLER": 4.114,
69  }),
70  "Y1S-off": ("Y1S", { # m(Y1S) - 60 MeV = 9.400 GeV
71  "energyHER": 6.223,
72  "energyLER": 3.556,
73  }),
74  "Y2S-off": ("Y2S", { # m(Y2S) - 60 MeV = 9.963 GeV
75  "energyHER": 6.596,
76  "energyLER": 3.769,
77  }),
78  "Y3S-off": ("Y3S", { # m(Y3S) - 60 MeV = 10.295 GeV
79  "energyHER": 6.816,
80  "energyLER": 3.895,
81  }),
82  "Y4S-off": ("Y4S", { # m(Y4S) - 60 MeV = 10.519 GeV
83  "energyHER": 6.964,
84  "energyLER": 3.979,
85  }),
86  "Y5S-off": ("Y5S", { # m(Y5S) - 60 MeV = 10.816 GeV
87  "energyHER": 7.160,
88  "energyLER": 4.092,
89  }),
90  "Y1D1": ("SuperKEKB", { # m(Y1D) = 10.152 GeV
91  "energyHER": 6.720,
92  "energyLER": 3.840,
93  }),
94  "Y2D1": ("SuperKEKB", { # m(Y2D) = 10.425 GeV
95  "energyHER": 6.901,
96  "energyLER": 3.944,
97  }),
98  "Y6S": ("SuperKEKB", { # m(Y6S) = 11.020 GeV
99  "energyHER": 7.295,
100  "energyLER": 4.169,
101  }),
102  "Yb10775": ("SuperKEKB", { # m(Yb10775) = 10.775 GeV
103  "energyHER": 7.133,
104  "energyLER": 4.076,
105  }),
106  "KEKB-Belle": (None, {
107  "energyHER": 7.998213,
108  "energyLER": 3.499218,
109  "angleXHER": 0.022,
110  }),
111  "LEP1": (None, {
112  "energyHER": 45.6,
113  "energyLER": 45.6,
114  }),
115  "DAPHNE": (None, {
116  "energyHER": 1.02 / 2.,
117  "energyLER": 1.02 / 2.,
118  }),
119 }
120 
121 
122 def __get_4vector(energy, angle):
123  """Calculate an 4vector for electron/positron from energy and angle"""
124  import ROOT
125  m = 0.511e-3
126  pz = (energy ** 2 - m ** 2) ** .5
127  v = ROOT.Math.PxPyPzEVector(0, 0, pz, energy)
128  if angle < 0:
129  angle = math.pi + angle
130  rotationY = ROOT.Math.RotationY(angle)
131  v = rotationY(v)
132  return v
133 
134 
135 def rot_matrix_y(angle):
136  """Return rotation matrix for a rotation of angle around the y-axis"""
137  c = math.cos(angle)
138  s = math.sin(angle)
139  return np.matrix([(c, 0, s), (0, 1, 0), (-s, 0, c)])
140 
141 
142 def cov_matrix(beamsize, angle):
143  """Return the covariance matrix for the given bunch size and angle with
144  respect to the z axis"""
145  if angle < 0:
146  angle = math.pi + angle
147  rot = rot_matrix_y(angle)
148  cov = np.matrix(np.diag(beamsize ** 2))
149  return rot * cov * rot.T
150 
151 
152 def calculate_beamspot(pos_her, pos_ler, size_her, size_ler, angle_her, angle_ler):
153  """
154  Function to calculate position and covariance matrix of the beamspot.
155 
156  This function calculates the mean position and the covariance matrix for the
157  beamspot from the bunch positions, bunch sizes and the angles of the beam
158  with respect to the z-axis. It assumes normal distributed bunch densities and
159  multiplies the PDFs for the two bunches to get a resulting beamspot
160 
161  negative angles are interpreted as "math.pi - abs(angle)"
162 
163  Args:
164  pos_her (list): mean position of the HER bunch (in cm)
165  pos_ler (list): mean position of the LER bunch (in cm)
166  size_her (list): HER bunch size in cm
167  size_ler (list): LER bunch size in cm
168  angle_her (float): angle between HER direction and z-axis
169  angle_her (float): angle between LER direction and z-axis
170 
171  Returns:
172  pos (array): array (x, y, z) with the spot center in cm
173  cov (matrix): covariance matrix of the beam spot
174  """
175 
176  cov_her = cov_matrix(size_her, angle_her)
177  cov_ler = cov_matrix(size_ler, angle_ler)
178  # calculate the inverse of the sum of the two covariance matrices
179  cov_sum_i = (cov_her + cov_ler).I
180  # and use it to calculate the product of the two gaussian distributions,
181  # covariance and mean
182  cov_spot = cov_her * cov_sum_i * cov_ler
183  beampos_spot = np.array(
184  cov_ler * cov_sum_i * np.matrix(pos_her).T +
185  cov_her * cov_sum_i * np.matrix(pos_ler).T
186  ).T[0]
187  return beampos_spot, cov_spot
188 
189 
190 def add_beamparameters(path, name, E_cms=None, **argk):
191  """Add BeamParameter module to a given path
192 
193  Args:
194  path (basf2.Path instance): path to add the module to
195  name (str): name of the beamparameter settings to use
196  E_cms (float): center of mass energy. If not None the beamenergies will
197  be scaled accordingly to achieve E_cms
198 
199  Additional keyword arguments will be passed directly to the module as parameters.
200  """
201  # calculate all the parameter values from the preset
202  values = calculate_beamparameters(name, E_cms)
203  # add a BeamParametes module to the path
204  module = path.add_module("BeamParameters")
205  module.set_name("BeamParameters:%s" % name)
206  # finally, set all parameters and return the module
207  module.param(values)
208  # and override parameters with any additional keyword arguments
209  module.param(argk)
210  return module
211 
212 
213 def calculate_beamparameters(name, E_cms=None):
214  """Get/calculate the necessary beamparameter values from the given preset name,
215  optionally scale it to achieve the given cms energy
216 
217  Args:
218  name (str): name of the beamparameter settings to use
219  E_cms (float): center of mass energy. If not None the beamenergies will
220  be scaled accordingly to achieve E_cms
221  """
222  if name not in beamparameter_presets:
223  B2FATAL("Unknown beamparameter preset: '%s', use one of %s" %
224  (name, ", ".join(sorted(beamparameter_presets.keys()))))
225  return None
226 
227  # copy the defaults
228  values = beamparameter_defaults.copy()
229 
230  # collect all the values from the preset by also looping over the
231  preset = {}
232 
233  def collect_values(preset, name):
234  parent, values = beamparameter_presets[name]
235  if parent is not None:
236  collect_values(preset, parent)
237  preset.update(values)
238 
239  collect_values(preset, name)
240 
241  # and update with the preset
242  for k in values.keys():
243  if k in preset:
244  if isinstance(preset[k], np.ndarray):
245  values[k] = list(preset[k].flat)
246  else:
247  values[k] = preset[k]
248 
249  # if we have bunch sizes we want to compute the vertex covariance from them
250  if "bunchHER" in preset and "bunchLER" in preset:
251  her_size = preset["bunchHER"]
252  ler_size = preset["bunchLER"]
253  her_angle = values["angleXHER"]
254  ler_angle = values["angleXLER"]
255  _, cov = calculate_beamspot(
256  [0, 0, 0], [0, 0, 0], her_size, ler_size, her_angle, ler_angle
257  )
258  values["covVertex"] = list(cov.flat)
259 
260  if E_cms is not None:
261  ler = __get_4vector(values["energyLER"], values["angleXLER"])
262  her = __get_4vector(values["energyHER"], values["angleXHER"])
263  mass = (her + ler).M()
264  scale = E_cms / mass
265  B2INFO("Scaling beam energies by %g to obtain E_cms = %g GeV" % (scale, E_cms))
266  values["energyHER"] *= scale
267  values["energyLER"] *= scale
268 
269  return values
270 
271 
272 if __name__ == "__main__":
273  # if called directly we will
274  # 1. calculate the beamspot for the SuperKEKB preset and display it if possible
275  # 2. calculate beam energies for Y1S - Y5S and offresonance
276  np.set_printoptions(precision=3)
277 
278  # calculate beamspot for SuperKEKB for testing
279 
280  values = beamparameter_presets["SuperKEKB"][1]
281  beampos_spot, cov_spot = calculate_beamspot(
282  [0, 0, 0], [0, 0, 0],
283  values["bunchHER"], values["bunchLER"],
284  values["angleXHER"], values["angleXLER"],
285  )
286  print("Beamspot position:")
287  print(beampos_spot)
288  print("Beamspot covariance:")
289  print(cov_spot)
290  print("Beamspot dimensions (in mm, axes in arbitary order):")
291  print(np.linalg.eig(cov_spot)[0] ** .5 * 10)
292 
293  # see if we can plot it
294  try:
295  from matplotlib import pyplot as pl
296 
297  cov_her = cov_matrix(values["bunchHER"], values["angleXHER"])
298 
299  cov_ler = cov_matrix(values["bunchLER"], values["angleXLER"])
300 
301  points_her = np.random.multivariate_normal([0, 0, 0], cov_her, 1000)
302 
303  points_ler = np.random.multivariate_normal([0, 0, 0], cov_ler, 1000)
304 
305  points_spot = np.random.multivariate_normal(beampos_spot, cov_spot, 1000)
306  pl.scatter(points_her[:, 2], points_her[:, 0], s=5, marker=".", c="b", edgecolors="None", label="HER")
307  pl.scatter(points_ler[:, 2], points_ler[:, 0], s=5, marker=".", c="r", edgecolors="None", label="LER")
308  pl.scatter(points_spot[:, 2], points_spot[:, 0], s=5, marker=".", c="g", edgecolors="None", label="spot")
309  pl.legend()
310  pl.xlabel("z / cm")
311  pl.ylabel("x / cm")
312  pl.tight_layout()
313  pl.show()
314  except ImportError:
315  pass
316 
317 
318  eher = values["energyHER"]
319 
320  eler = values["energyLER"]
321 
322  aher = values["angleXHER"]
323 
324  aler = values["angleXLER"]
325 
326  # calculate beam energies for Y1S - Y4S by scaling them from the nominal
327  # beam energies for the SuperKEKB preset
328 
329  ler = __get_4vector(eher, aher)
330 
331  her = __get_4vector(eler, aler)
332 
333  mass = (her + ler).M()
334  print("Nominal CMS Energy: ", mass)
335 
336 
337  targets = {
338  "Y1S": 9460.30e-3,
339  "Y2S": 10023.26e-3,
340  "Y3S": 10355.2e-3,
341  "Y4S": 10579.4e-3,
342  "Y5S": 10876e-3
343  }
344  for name, energy in sorted(targets.items()):
345 
346  scale = energy / mass
347  print("""\
348  "%s": ("SuperKEKB", { # m(%s) = %.3f GeV
349  "energyHER": %.3f,
350  "energyLER": %.3f,
351  }),""" % (name, name, energy, eher * scale, eler * scale))
352 
353  for name, energy in sorted(targets.items()):
354 
355  scale = (energy - 60e-3) / mass
356  print("""\
357  "%s-off": ("%s", { # m(%s) - 60 MeV = %.3f GeV
358  "energyHER": %.3f,
359  "energyLER": %.3f,
360  }),""" % (name, name, name, energy - 60e-3, eher * scale, eler * scale))