13 Scripts and functions to set the BeamParameters from known configurations
16 from basf2
import B2FATAL, B2INFO
21 beamparameter_defaults = {
40 beamparameter_presets = {
52 "covHER": np.array([0.00513]) ** 2,
53 "covLER": np.array([0.002375]) ** 2,
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,
58 "Y1S": (
"SuperKEKB", {
62 "Y2S": (
"SuperKEKB", {
66 "Y3S": (
"SuperKEKB", {
70 "Y4S": (
"SuperKEKB", {
74 "Y5S": (
"SuperKEKB", {
98 "Y1D1": (
"SuperKEKB", {
102 "Y2D1": (
"SuperKEKB", {
106 "Y6S": (
"SuperKEKB", {
110 "Yb10775": (
"SuperKEKB", {
114 "KEKB-Belle": (
None, {
115 "energyHER": 7.998213,
116 "energyLER": 3.499218,
124 "energyHER": 1.02 / 2.,
125 "energyLER": 1.02 / 2.,
130 def __get_4vector(energy, angle):
131 """Calculate an 4vector for electron/positron from energy and angle"""
134 pz = (energy ** 2 - m ** 2) ** .5
135 v = ROOT.Math.PxPyPzEVector(0, 0, pz, energy)
137 angle = math.pi + angle
138 rotationY = ROOT.Math.RotationY(angle)
143 def rot_matrix_y(angle):
144 """Return rotation matrix for a rotation of angle around the y-axis"""
147 return np.matrix([(c, 0, s), (0, 1, 0), (-s, 0, c)])
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"""
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
160 def calculate_beamspot(pos_her, pos_ler, size_her, size_ler, angle_her, angle_ler):
162 Function to calculate position and covariance matrix of the beamspot.
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
169 negative angles are interpreted as "math.pi - abs(angle)"
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
180 pos (array): array (x, y, z) with the spot center in cm
181 cov (matrix): covariance matrix of the beam spot
184 cov_her = cov_matrix(size_her, angle_her)
185 cov_ler = cov_matrix(size_ler, angle_ler)
187 cov_sum_i = (cov_her + cov_ler).I
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
195 return beampos_spot, cov_spot
198 def add_beamparameters(path, name, E_cms=None, **argk):
199 """Add BeamParameter module to a given path
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
207 Additional keyword arguments will be passed directly to the module as parameters.
210 values = calculate_beamparameters(name, E_cms)
212 module = path.add_module(
"BeamParameters")
213 module.set_name(
"BeamParameters:%s" % name)
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
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
230 if name
not in beamparameter_presets:
231 B2FATAL(
"Unknown beamparameter preset: '%s', use one of %s" %
232 (name,
", ".join(sorted(beamparameter_presets.keys()))))
236 values = beamparameter_defaults.copy()
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)
247 collect_values(preset, name)
250 for k
in values.keys():
252 if isinstance(preset[k], np.ndarray):
253 values[k] = list(preset[k].flat)
255 values[k] = preset[k]
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
266 values[
"covVertex"] = list(cov.flat)
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()
273 B2INFO(
"Scaling beam energies by %g to obtain E_cms = %g GeV" % (scale, E_cms))
274 values[
"energyHER"] *= scale
275 values[
"energyLER"] *= scale
280 def _get_collisions_invariant_mass(experiment, run, verbose=False):
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.
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
291 from ROOT
import Belle2
as B2
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)
299 beams = B2.PyDBObj(
'BeamParameters')
300 if not beams.isValid():
301 B2FATAL(
'BeamParameters is not valid')
302 collisions_invariant_mass = beams.getMass()
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()}')
312 B2.DataStore.Instance().reset()
314 return collisions_invariant_mass
317 def get_collisions_invariant_mass(experiment, run, verbose=False):
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.
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
329 def _run_function_in_process(func, *args, **kwargs):
331 Small helper function for running another function in a process and thus avoiding undesired side effects.
333 from multiprocessing
import Pool
335 with Pool(processes=1)
as pool:
337 result = pool.apply_async(func, args, kwargs)
339 output = result.get()
343 return _run_function_in_process(
344 func=_get_collisions_invariant_mass, experiment=experiment, run=run, verbose=verbose
348 if __name__ ==
"__main__":
352 np.set_printoptions(precision=3)
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"],
362 print(
"Beamspot position:")
364 print(
"Beamspot covariance:")
366 print(
"Beamspot dimensions (in mm, axes in arbitary order):")
367 print(np.linalg.eig(cov_spot)[0] ** .5 * 10)
371 from matplotlib
import pyplot
as pl
373 cov_her = cov_matrix(values[
"bunchHER"], values[
"angleXHER"])
375 cov_ler = cov_matrix(values[
"bunchLER"], values[
"angleXLER"])
377 points_her = np.random.multivariate_normal([0, 0, 0], cov_her, 1000)
379 points_ler = np.random.multivariate_normal([0, 0, 0], cov_ler, 1000)
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")
394 eher = values[
"energyHER"]
396 eler = values[
"energyLER"]
398 aher = values[
"angleXHER"]
400 aler = values[
"angleXLER"]
405 ler = __get_4vector(eher, aher)
407 her = __get_4vector(eler, aler)
409 mass = (her + ler).M()
410 print(
"Nominal CMS Energy: ", mass)
420 for name, energy
in sorted(targets.items()):
422 scale = energy / mass
424 "%s": ("SuperKEKB", { # m(%s) = %.3f GeV
427 }),""" % (name, name, energy, eher * scale, eler * scale))
429 for name, energy
in sorted(targets.items()):
431 scale = (energy - 60e-3) / mass
433 "%s-off": ("%s", { # m(%s) - 60 MeV = %.3f GeV
436 }),""" % (name, name, name, energy - 60e-3, eher * scale, eler * scale))