Belle II Software  release-06-00-14
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.TLorentzVector(0, 0, pz, energy)
128  if angle < 0:
129  angle = math.pi + angle
130  v.RotateY(angle)
131  return v
132 
133 
134 def rot_matrix_y(angle):
135  """Return rotation matrix for a rotation of angle around the y-axis"""
136  c = math.cos(angle)
137  s = math.sin(angle)
138  return np.matrix([(c, 0, s), (0, 1, 0), (-s, 0, c)])
139 
140 
141 def cov_matrix(beamsize, angle):
142  """Return the covariance matrix for the given bunch size and angle with
143  respect to the z axis"""
144  if angle < 0:
145  angle = math.pi + angle
146  rot = rot_matrix_y(angle)
147  cov = np.matrix(np.diag(beamsize ** 2))
148  return rot * cov * rot.T
149 
150 
151 def calculate_beamspot(pos_her, pos_ler, size_her, size_ler, angle_her, angle_ler):
152  """
153  Function to calculate position and covariance matrix of the beamspot.
154 
155  This function calculates the mean position and the covariance matrix for the
156  beamspot from the bunch positions, bunch sizes and the angles of the beam
157  with respect to the z-axis. It assumes normal distributed bunch densities and
158  multiplies the PDFs for the two bunches to get a resulting beamspot
159 
160  negative angles are interpreted as "math.pi - abs(angle)"
161 
162  Args:
163  pos_her (list): mean position of the HER bunch (in cm)
164  pos_ler (list): mean position of the LER bunch (in cm)
165  size_her (list): HER bunch size in cm
166  size_ler (list): LER bunch size in cm
167  angle_her (float): angle between HER direction and z-axis
168  angle_her (float): angle between LER direction and z-axis
169 
170  Returns:
171  pos (array): array (x, y, z) with the spot center in cm
172  cov (matrix): covariance matrix of the beam spot
173  """
174 
175  cov_her = cov_matrix(size_her, angle_her)
176  cov_ler = cov_matrix(size_ler, angle_ler)
177  # calculate the inverse of the sum of the two covariance matrices
178  cov_sum_i = (cov_her + cov_ler).I
179  # and use it to calculate the product of the two gaussian distributions,
180  # covariance and mean
181  cov_spot = cov_her * cov_sum_i * cov_ler
182  beampos_spot = np.array(
183  cov_ler * cov_sum_i * np.matrix(pos_her).T +
184  cov_her * cov_sum_i * np.matrix(pos_ler).T
185  ).T[0]
186  return beampos_spot, cov_spot
187 
188 
189 def add_beamparameters(path, name, E_cms=None, **argk):
190  """Add BeamParameter module to a given path
191 
192  Args:
193  path (basf2.Path instance): path to add the module to
194  name (str): name of the beamparameter settings to use
195  E_cms (float): center of mass energy. If not None the beamenergies will
196  be scaled accordingly to achieve E_cms
197 
198  Additional keyword arguments will be passed directly to the module as parameters.
199  """
200  # calculate all the parameter values from the preset
201  values = calculate_beamparameters(name, E_cms)
202  # add a BeamParametes module to the path
203  module = path.add_module("BeamParameters")
204  module.set_name("BeamParameters:%s" % name)
205  # finally, set all parameters and return the module
206  module.param(values)
207  # and override parameters with any additional keyword arguments
208  module.param(argk)
209  return module
210 
211 
212 def calculate_beamparameters(name, E_cms=None):
213  """Get/calculate the necessary beamparameter values from the given preset name,
214  optionally scale it to achieve the given cms energy
215 
216  Args:
217  name (str): name of the beamparameter settings to use
218  E_cms (float): center of mass energy. If not None the beamenergies will
219  be scaled accordingly to achieve E_cms
220  """
221  if name not in beamparameter_presets:
222  B2FATAL("Unknown beamparameter preset: '%s', use one of %s" %
223  (name, ", ".join(sorted(beamparameter_presets.keys()))))
224  return None
225 
226  # copy the defaults
227  values = beamparameter_defaults.copy()
228 
229  # collect all the values from the preset by also looping over the
230  preset = {}
231 
232  def collect_values(preset, name):
233  parent, values = beamparameter_presets[name]
234  if parent is not None:
235  collect_values(preset, parent)
236  preset.update(values)
237 
238  collect_values(preset, name)
239 
240  # and update with the preset
241  for k in values.keys():
242  if k in preset:
243  if isinstance(preset[k], np.ndarray):
244  values[k] = list(preset[k].flat)
245  else:
246  values[k] = preset[k]
247 
248  # if we have bunch sizes we want to compute the vertex covariance from them
249  if "bunchHER" in preset and "bunchLER" in preset:
250  her_size = preset["bunchHER"]
251  ler_size = preset["bunchLER"]
252  her_angle = values["angleXHER"]
253  ler_angle = values["angleXLER"]
254  _, cov = calculate_beamspot(
255  [0, 0, 0], [0, 0, 0], her_size, ler_size, her_angle, ler_angle
256  )
257  values["covVertex"] = list(cov.flat)
258 
259  if E_cms is not None:
260  ler = __get_4vector(values["energyLER"], values["angleXLER"])
261  her = __get_4vector(values["energyHER"], values["angleXHER"])
262  mass = (her + ler).M()
263  scale = E_cms / mass
264  B2INFO("Scaling beam energies by %g to obtain E_cms = %g GeV" % (scale, E_cms))
265  values["energyHER"] *= scale
266  values["energyLER"] *= scale
267 
268  return values
269 
270 
271 if __name__ == "__main__":
272  # if called directly we will
273  # 1. calculate the beamspot for the SuperKEKB preset and display it if possible
274  # 2. calculate beam energies for Y1S - Y5S and offresonance
275  np.set_printoptions(precision=3)
276 
277  # calculate beamspot for SuperKEKB for testing
278 
279  values = beamparameter_presets["SuperKEKB"][1]
280  beampos_spot, cov_spot = calculate_beamspot(
281  [0, 0, 0], [0, 0, 0],
282  values["bunchHER"], values["bunchLER"],
283  values["angleXHER"], values["angleXLER"],
284  )
285  print("Beamspot position:")
286  print(beampos_spot)
287  print("Beamspot covariance:")
288  print(cov_spot)
289  print("Beamspot dimensions (in mm, axes in arbitary order):")
290  print(np.linalg.eig(cov_spot)[0] ** .5 * 10)
291 
292  # see if we can plot it
293  try:
294  from matplotlib import pyplot as pl
295 
296  cov_her = cov_matrix(values["bunchHER"], values["angleXHER"])
297 
298  cov_ler = cov_matrix(values["bunchLER"], values["angleXLER"])
299 
300  points_her = np.random.multivariate_normal([0, 0, 0], cov_her, 1000)
301 
302  points_ler = np.random.multivariate_normal([0, 0, 0], cov_ler, 1000)
303 
304  points_spot = np.random.multivariate_normal(beampos_spot, cov_spot, 1000)
305  pl.scatter(points_her[:, 2], points_her[:, 0], s=5, marker=".", c="b", edgecolors="None", label="HER")
306  pl.scatter(points_ler[:, 2], points_ler[:, 0], s=5, marker=".", c="r", edgecolors="None", label="LER")
307  pl.scatter(points_spot[:, 2], points_spot[:, 0], s=5, marker=".", c="g", edgecolors="None", label="spot")
308  pl.legend()
309  pl.xlabel("z / cm")
310  pl.ylabel("x / cm")
311  pl.tight_layout()
312  pl.show()
313  except ImportError:
314  pass
315 
316 
317  eher = values["energyHER"]
318 
319  eler = values["energyLER"]
320 
321  aher = values["angleXHER"]
322 
323  aler = values["angleXLER"]
324 
325  # calculate beam energies for Y1S - Y4S by scaling them from the nominal
326  # beam energies for the SuperKEKB preset
327 
328  ler = __get_4vector(eher, aher)
329 
330  her = __get_4vector(eler, aler)
331 
332  mass = (her + ler).M()
333  print("Nominal CMS Energy: ", mass)
334 
335 
336  targets = {
337  "Y1S": 9460.30e-3,
338  "Y2S": 10023.26e-3,
339  "Y3S": 10355.2e-3,
340  "Y4S": 10579.4e-3,
341  "Y5S": 10876e-3
342  }
343  for name, energy in sorted(targets.items()):
344 
345  scale = energy / mass
346  print("""\
347  "%s": ("SuperKEKB", { # m(%s) = %.3f GeV
348  "energyHER": %.3f,
349  "energyLER": %.3f,
350  }),""" % (name, name, energy, eher * scale, eler * scale))
351 
352  for name, energy in sorted(targets.items()):
353 
354  scale = (energy - 60e-3) / mass
355  print("""\
356  "%s-off": ("%s", { # m(%s) - 60 MeV = %.3f GeV
357  "energyHER": %.3f,
358  "energyLER": %.3f,
359  }),""" % (name, name, name, energy - 60e-3, eher * scale, eler * scale))