Optics
Beam Optics
Module implementing various functionality for simple beam parameter calculations.
- class pyhdtoolkit.optics.beam.Beam(energy: float, gemitt: float, m0: float = 938.27208943)[source]
- Added in version 0.6.0. - Class to represent most useful particle beam attributes for - MAD-Xsimulations.- property beta_rel: float
- Relativistic beta. 
 - property brho: float
- Beam rigidity in [T/m]. 
 - eta(alpha_p: float) float[source]
- Returns the slip factor parameter \(\eta\). - Note - We use the convention according to which \(\eta = 0\) at transition energy and \(\eta < 0\) above transition. 
 - property gamma_rel: float
- Relativistic gamma. 
 - static gamma_transition(alpha_p: float) float[source]
- Returns the relativistic \(\gamma\) corresponding to the transition energy. 
 - property nemitt: float
- Normalized emittance in [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. 
 - property rms_emittance: float
- Rms emittance in [m]. 
 
- pyhdtoolkit.optics.beam.compute_beam_parameters(pc_gev: float, nemitt_x: float, nemitt_y: float, deltap_p: float) BeamParameters[source]
- Added 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 gallery. - Parameters:
- Returns:
- BeamParameters-- A- BeamParametersobject 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 
 
 
Resonance Driving Terms Utilities
Module implementing utilities for the handling of resonance driving terms.
- pyhdtoolkit.optics.rdt.determine_rdt_line(rdt: int | str, plane: str) tuple[int, int, int][source]
- Added in version 1.5.0. - Find the given line to look for in the spectral analysis of the given plane that corresponds to the given RDT. - Parameters:
- rdt ( - int | str) -- The RDT to look for.
- plane ( - str) -- The plane to look for the RDT in.
 
- Returns:
- tuple[int,- int,- int]-- A tuple of three integers representing the line to look for in the spectral analysis. For instance, f1001 corresponds to the line (0, 1, 0) in the X plane which means the line located at 1 * Qy = Qy.
 - Examples - determine_rdt_line(1001, "X") # Output: (0, 1, 0) # Line at 1 * Qy in the X spectrum. - determine_rdt_line("2002", "Y") # Output: (-2, 3, 0) # Line at 3 * Qy - 2 * Qx in the Y spectrum. 
- pyhdtoolkit.optics.rdt.rdt_to_order_and_type(rdt: int | str) str[source]
- Added in version 1.5.0. - Decompose the input RDT into its four various components and return the type of RDT (normal or skew) and its order. - Parameters:
- rdt ( - int | str) -- The RDT to decompose.
- Returns:
- str-- A string with the type and magnet order of the provided RDT.
 - Examples - rdt_to_order_and_type(1001) # Output: 'skew_quadrupole' - rdt_to_order_and_type("2002") # Output: 'normal_octupole' 
Ripken Parameters
Module implementing various calculations based on the Willeke and Ripken [WR89] optics parameters.
- pyhdtoolkit.optics.ripken.lebedev_beam_size(beta1_: float | ndarray, beta2_: float | ndarray, gemitt_x: float, gemitt_y: float) float | ndarray[source]
- Added 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]. - Hint - 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:
- Returns:
- float | numpy.ndarray-- The beam size (horizontal or vertical) according to Lebedev& Bogacz, as \(\sqrt{\epsilon_x * \beta_{1,\_}^2 + \epsilon_y * \beta_{2,\_}^2}\).
 - Example - gemitt_x = madx.globals["geometric_emit_x"] gemitt_y = madx.globals["geometric_emit_y"] twiss_tfs = madx.twiss(ripken=True).dframe() horizontal_size = lebedev_beam_size( twiss_tfs.beta11, twiss_tfs.beta21, gemitt_x, gemitt_y ) vertical_size = lebedev_beam_size( twiss_tfs.beta12, twiss_tfs.beta22, gemitt_x, gemitt_y ) 
Twiss Optics
Module implementing various calculations based on the TWISS optics parameters.
- pyhdtoolkit.optics.twiss.courant_snyder_transform(u_vector: ndarray, alpha: float, beta: float) ndarray[source]
- Added 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 ( - numpy.ndarray) -- Two-dimentional array of the phase-space (spatial and momentum) coordinates, either horizontal or vertical.
- alpha ( - float) -- Alpha twiss parameter in the appropriate plane.
- beta ( - float) -- Beta twiss parameter in the appropriate plane.
 
- Returns:
- numpy.ndarray-- 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)