Optics

Beam Optics

Module implementing various functionality for simple beam parameter calculations.

class pyhdtoolkit.optics.beam.Beam(energy: float, emittance: float, m0: float = 938.27208816)[source]

New in version 0.6.0.

Class to represent most useful particle beam attributes for MAD-X simulations.

__init__(energy: float, emittance: float, m0: float = 938.27208816) None[source]
Parameters
  • energy (float) -- energy of the particles in your beam, in [GeV].

  • emittance (float) -- beam emittance, in [m].

  • m0 (float) -- rest mass of the beam’s particles in [MeV]. Defaults to that of a proton.

property beta_rel: float

Relativistic beta.

property brho: float

Beam rigidity [T/m].

eta(alpha_p: float) float[source]

Returns the slip factor parameter \(\eta\).

Note

\(\eta = 0\) at transition energy (\(\eta < 0\) above transition).

Parameters

alpha_p (float) -- momentum compaction factor.

Returns

The slip factor.

property gamma_rel: float

Relativistic gamma.

static gamma_transition(alpha_p: float) float[source]

Returns the relativistic \(\gamma\) corresponding to the transition energy.

Parameters

alpha_p (float) -- momentum compaction factor.

Returns

The relativistic \(\gamma\) value at the transition energy.

property normalized_emittance: float

Normalized emittance [m].

revolution_frequency(circumference: float = 26658.8832, speed: float = 299792458.0) float[source]

Returns the revolution frequency of the beam’s particles around the accelerator.

Parameters
  • circumference (float) -- the machine circumference in [m]. Defaults to that of the LHC.

  • speed (float) -- the particles’ speed in the machine, in [m/s]. Defaults to c.

Returns

The revolution frequency, in [turns/s].

property rms_emittance: float

Rms emittance [m].

pyhdtoolkit.optics.beam.compute_beam_parameters(pc_gev: float, en_x_m: float, en_y_m: float, deltap_p: float) pyhdtoolkit.models.beam.BeamParameters[source]

New in version 0.12.0.

Calculates beam parameters from provided values, for proton particles. One can find an example use of this function in the beam enveloppe example gallery.

Parameters
  • pc_gev (float) -- particle momentum, in [GeV].

  • en_x_m (float) -- horizontal emittance, in [m].

  • en_y_m (float) -- vertical emittance, in [m].

  • deltap_p (float) -- momentum deviation.

Returns

A BeamParameters object with the calculated values.

Example

params = compute_beam_parameters(1.9, 5e-6, 5e-6, 2e-3)
print(params)
# Beam Parameters for particle of charge 1
# Beam momentum = 1.900 GeV/c
# Normalized x-emittance = 5.000 mm mrad
# Normalized y-emittance = 5.000 mm mrad
# Momentum deviation deltap/p = 0.002
# -> Beam total energy = 2.119 GeV
# -> Beam kinetic energy = 1.181 GeV
# -> Beam rigidity = 6.333 Tm
# -> Relativistic beta = 0.89663
# -> Relativistic gamma = 2.258
# -> Geometrical x emittance = 2.469 mm mrad
# -> Geometrical y emittance = 2.469 mm mrad

Ripken Parameters

Module implementing various calculations based on the Willeke and Ripken [WR89] optics parameters.

pyhdtoolkit.optics.ripken.lebedev_beam_size(beta1_: Union[float, numpy.ndarray], beta2_: Union[float, numpy.ndarray], geom_emit_x: float, geom_emit_y: float) Union[float, numpy.ndarray][source]

New in version 0.8.2.

Calculate beam size according to the Lebedev-Bogacz formula, based on the Ripken-Mais Twiss parameters. The implementation is that of Eq. (A.3.1) in Lebedev and Bogacz [LB10].

Tip

For the calculations, use \(\beta_{11}\) and \(\beta_{21}\) for the vertical beam size, but use \(\beta_{12}\) and \(\beta_{22}\) for the horizontal one. See the example below.

Parameters
  • beta1 (Union[float, np.ndarray]) -- value(s) for the beta1x or beta1y Ripken parameter.

  • beta2 (Union[float, np.ndarray]) -- value(s) for the beta2x or beta2y Ripken parameter.

  • geom_emit_x (float) -- geometric emittance of the horizontal plane, in [m].

  • geom_emit_y (float) -- geometric emittante of the vertical plane, in [m].

Returns

The beam size (horizontal or vertical) according to Lebedev & Bogacz, as \(\sqrt{\epsilon_x * \beta_{1,\_}^2 + \epsilon_y * \beta_{2,\_}^2}\).

Example

geom_emit_x = madx.globals["geometric_emit_x"]
geom_emit_y = madx.globals["geometric_emit_y"]
twiss_tfs = madx.twiss(ripken=True).dframe()
horizontal_size = lebedev_beam_size(
    twiss_tfs.beta11, twiss_tfs.beta21, geom_emit_x, geom_emit_y
)
vertical_size = lebedev_beam_size(
    twiss_tfs.beta12, twiss_tfs.beta22, geom_emit_x, geom_emit_y
)

Twiss Optics

Module implementing various calculations based on the TWISS optics parameters.

pyhdtoolkit.optics.twiss.courant_snyder_transform(u_vector: numpy.ndarray, alpha: float, beta: float) numpy.ndarray[source]

New in version 0.5.0.

Perform the Courant-Snyder transform on regular (nonchaotic) phase-space coordinates. Specifically, if considering the horizontal plane and noting \(U = (x, px)\) the phase-space vector, it returns \(\bar{U} = (\bar{x}, \bar{px})\) according to the transform \(\bar{U} = P \cdot U\), where

\[\begin{split}P = \begin{pmatrix} \frac{1}{\sqrt{\beta_x}} & 0 \\ \frac{\alpha_x}{\sqrt{\beta_x}} & \sqrt{\beta_x} \\ \end{pmatrix}\end{split}\]
Parameters
  • u_vector (np.ndarray) -- two-dimentional array of phase-space (spatial and momenta) coordinates, either horizontal or vertical.

  • alpha (float) -- alpha twiss parameter in the appropriate plane.

  • beta (float) -- beta twiss parameter in the appropriate plane.

Returns

The normalized phase-space coordinates from the Courant-Snyder transform.

Example

alfx = madx.table.twiss.alfx[0]
betx = madx.table.twiss.betx[0]
u = np.array([x_coords, px_coord])
u_bar = courant_snyder_transform(u, alfx, betx)