Parameters: Python

This documents on how to use ImpactX as a Python script (python3 run_script.py).

Tip

If you enjoy AI/LLM/agentic workflows, see our AI (LLM)-Assisted Input File Design section, too.

Collective Effects & Overall Simulation Parameters

class impactx.ImpactX

This is the central simulation class.

property particle_shape

Control the particle B-spline order.

The order of the shape factors (splines) for the macro-particles along all spatial directions: 1 for linear, 2 for quadratic, 3 for cubic. Low-order shape factors result in faster simulations, but may lead to more noisy results. High-order shape factors are computationally more expensive, but may increase the overall accuracy of the results. For production runs it is generally safer to use high-order shape factors, such as cubic order.

Parameters:

order (int) – B-spline order 1, 2, or 3

property n_cell

The number of grid points along each direction on the coarsest level.

property max_level

The maximum mesh-refinement level for the simulation.

property finest_level

The currently finest level of mesh-refinement used. This is always less or equal to max_level.

property domain

The physical extent of the full simulation domain, relative to the reference particle position, in meters. When set, turns dynamic_size to False.

Note: particles that move outside the simulation domain are removed.

property prob_relative

This is a list with amr.max_level + 1 entries.

By default, we dynamically extract the minimum and maximum of the particle positions in the beam. The field mesh spans, per direction, multiple times the maximum physical extent of beam particles, as given by this factor. The beam minimum and maximum extent are symmetrically padded by the mesh. For instance, 1.2 means the mesh will span 10% above and 10% below the beam; 1.0 means the beam is exactly covered with the mesh.

Default: [3.0 1.0 1.0 ...]. When set, turns dynamic_size to True.

property dynamic_size

Use dynamic (True) resizing of the field mesh or static sizing (False).

property space_charge

The physical model of space charge used.

Options:

  • False (default): space charge effects are not calculated.

  • "2D": Space charge forces are computed in the plane (x,y) transverse to the reference particle velocity, assuming the beam is long and unbunched.

  • "2p5D": Space charge forces are computed in the plane (x,y) transverse to the reference particle velocity, while the transverse space charge kicks are weighted by the longitudinal line density determined by charge deposition (2.5D model). Longitudinal space charge kicks are determined by the derivative of the line charge density.

  • "3D": Space charge forces are computed in three dimensions, assuming the beam is bunched.

    When running in envelope mode (when algo.track = "envelope"), this model currently assumes that <xy> = <yt> = <tx> = 0.

  • "Gauss3D": Calculate 3D space charge forces as if the beam was a Gaussian distribution.

  • "Gauss2p5D": Calculate 2.5D space charge forces as if the beam was a transverse Gaussian distribution.

    These models are supported only in particle tracking mode (when algo.track = "particles"). Ref.: J. Qiang, “Two-and-a-half dimensional symplectic space-charge solver”, LBNL Report Number: LBNL-2001674 (2025). (This reference describes both 3D and 2.5D models.)

property space_charge_gauss_nint

Number of steps for computing the integrals (default: 101).

property space_charge_gauss_taylor_delta

Initial integral region to avoid integrand divergence at 0 (default: 0.01).

property space_charge_gauss_charge_z_bins

Number of bins for longitudinal charge density deposition (default: 129). Used by the Gauss2p5D space charge model.

property space_charge_num_longitudinal_bins

Number of bins for longitudinal charge density deposition (default: 100). Used by the 2p5D space charge model.

property space_charge_apply_longitudinal_kick

Enable or disable the longitudinal space charge kick (default: True).

property poisson_solver

The numerical solver to solve the Poisson equation when calculating space charge effects. Either "fft" (default) or "multigrid".

Currently, the multigrid solver supports only 3D space charge. The fft solver supports 2D, 2.5D or 3D space charge.

  • fft: Poisson’s equation is solved using an Integrated Green Function method (which requires FFT calculations). See these references for more details Qiang et al. (2006) (+ Erratum). This requires the compilation flag -DImpactX_FFT=ON. If mesh refinement (MR) is enabled, this FFT solver is used only on the coarsest level and a multi-grid solver is used on refined levels. The boundary conditions are assumed to be open.

  • multigrid: Poisson’s equation is solved using an iterative multigrid (MLMG) solver. See the AMReX documentation for details of the MLMG solver.

property mlmg_relative_tolerance

Default: 1.e-7

The relative precision with which the electrostatic space-charge fields should be calculated. More specifically, the space-charge fields are computed with an iterative Multi-Level Multi-Grid (MLMG) solver. This solver can fail to reach the default precision within a reasonable time.

property mlmg_absolute_tolerance

Default: 0, which means: ignored

The absolute tolerance with which the space-charge fields should be calculated in units of \(V/m^2\). More specifically, the acceptable residual with which the solution can be considered converged. In general this should be left as the default, but in cases where the simulation state changes very little between steps it can occur that the initial guess for the MLMG solver is so close to the converged value that it fails to improve that solution sufficiently to reach the mlmg_relative_tolerance value.

property mlmg_max_iters

Default: 100

Maximum number of iterations used for MLMG solver for space-charge fields calculation. In case if MLMG converges but fails to reach the desired self_fields_required_precision, this parameter may be increased.

property mlmg_verbosity

Default: 1

The verbosity used for MLMG solver for space-charge fields calculation. Currently MLMG solver looks for verbosity levels from 0-5. A higher number results in more verbose output.

property csr

Enable (True) or disable (False) Coherent Synchrotron Radiation (CSR) calculations (default: False).

Whether to calculate Coherent Synchrotron Radiation (CSR) effects (default: disabled).

Currently, this is the 1D ultrarelativistic steady-state wakefield model (eq. 19 of E. L. Saldin et al, NIMA 398, p. 373-394 (1997), DOI:10.1016/S0168-9002(97)00822-X).

Note

CSR effects are only calculated for lattice elements that include bending, such as Sbend, ExactSbend and CFbend.

CSR effects require the compilation flag -DImpactX_FFT=ON.

property csr_bins

The number of bins along the longitudinal direction used for the CSR calculations (default: 150).

property isr

Enable (True) or disable (False) Incoherent Synchrotron Radiation (ISR) calculations (default: False).

Whether to calculate Incoherent Synchrotron Radiation (ISR) effects (default: disabled).

ISR effects are included in the simulation for bend lattice elements such as Sbend and CFbend, and are especially important for electron or positron bunches at high energy. The effects of ISR include radiation reaction due to the stochastic emission of synchrotron radiation, resulting in mean energy loss and quantum excitation of the bunch. The model is based on:

F. Niel et al., Phys. Rev. E 97, 043209 (2018), DOI:10.1103/PhysRevE.97.043209

However, a Taylor expansion is used to evaluate the dependence on the quantum parameter \(\chi\). When algo.isr_order = 1, the model is equivalent to that described in:

J. M. Jowett, “Introductory Statistical Mechanics for Electron Storage Rings”, AIP Conf. Proc. 153, 864-970 (1987), DOI:10.1063/1.36374

Note

ISR effects are only calculated for lattice elements that include bending, such as Sbend, ExactSbend and CFbend.

property isr_order

The number of terms retained in the Taylor series for the functions \(g(\chi)\) and \(h(\chi)\) appearing in Niel et al, equations (25) and (41) describing quantum effects.

property isr_on_ref_part

Flag specifying whether ISR is to be applied to the reference particle. When sim.isr_on_ref_part = False, the reference particle does not lose energy due to radiation, and the mean energy of the beam particles will decrease. This option is natural if the lattice optics, magnet settings, etc. are chosen without accounting for radiative energy loss. When sim.isr_on_ref_part = True, the reference particle does lose energy due to radiation, and little centroid evolution is expected in the beam particles. This option is natural if the lattice optics, magnet settings, etc. are chosen to account for radiative energy loss.

property particle_bc

The application of a non-trivial boundary condition for particles is currently supported only in the longitudinal direction (the local direction of motion as defined by the reference particle). This parameter sets the type of particle boundary condition to be applied to the longitudinal coordinate, based on the value of parameter bucket_length.

Options:

  • "open" (default): no action is applied at the boundary.

  • "periodic": each particle’s longitudinal coordinate is treated modulo the value bucket_length (phase wrapping).

  • "absorbing": a particle whose longitudinal coordinate falls outside the boundary is declared lost.

  • "reflecting": a particle whose longitudinal coordinate crosses the boundary is reflected about the boundary, with reversed longitudinal momentum.

    Note

    The implementation works through linear order in the phase space variables. If you have the need for a more precise implementation of reflecting boundaries, please open an issue.

property spin

Enable (True) or disable (False) particle spin tracking (default: False).

Whether to track particle spin. Currently, the implementation of spin tracking is a work in progress, and this feature is not yet supported.

property diagnostics

Enable (True) or disable (False) diagnostics generally (default: True). Disabling this is mostly used for benchmarking.

property slice_step_diagnostics

Enable (True) or disable (False) diagnostics every slice step in elements (default: True).

By default, diagnostics is performed at the beginning and end of the simulation. Enabling this flag will write diagnostics every step and slice step.

property diag_file_min_digits

The minimum number of digits (default: 6) used for the step number appended to the diagnostic file names.

property particle_lost_diagnostics_backend

Diagnostics for particles lost in apertures. See the BeamMonitor element for backend values.

property eigenemittances

Enable (True) or disable (False) output of eigenemittances at every slice step in elements (default: False).

If this flag is enabled, the 3 eigenemittances of the 6D beam distribution are computed and written as diagnostics. This flag is disabled by default to reduce computational cost.

init_grids()

Initialize AMReX blocks/grids for domain decomposition & space charge mesh.

This must come first, before particle beams and lattice elements are initialized.

add_particles(charge_C, distr, npart, spinv=None)

Particle tracking mode: Generate and add n particles to the particle container. Note: Set the reference particle properties (charge, mass, energy) first.

Will also resize the geometry based on the updated particle distribution’s extent and then redistribute particles in according AMReX grid boxes.

Parameters:
  • charge_C (float) – bunch charge (C)

  • distr – distribution function to draw from (object from impactx.distribution)

  • npart (int) – number of particles to draw

  • spinv (SpinvMF) – optional spin distribution

init_envelope(ref, distr, intensity=None)

Envelope tracking mode: Create a 6x6 covariance matrix from a distribution and then initialize the simulation for envelope tracking relative to a reference particle.

Parameters:
  • ref – the reference particle (object from impactx.RefPart)

  • distr – distribution function (object from impactx.distribution)

  • intensity (float) – the beam intensity, given as bunch charge (C) for 3D or beam current (A) for 2D space charge

particle_container()

Access the beam particle container (impactx.ParticleContainer). Deprecated: use beam.

property beam

Access the beam particle container (impactx.ParticleContainer).

property lattice

Access the elements in the accelerator lattice. See impactx.elements for lattice elements.

property periods

The number of periods to repeat the lattice.

property omp_threads

Controls the number of OpenMP threads to use (ImpactX default: “nosmt”).

See the detailed AMReX docs for details in the accepted values.

property abort_on_warning_threshold

(optional) Set to “low”, “medium” or “high”. Cause the code to abort if a warning is raised that exceeds the warning threshold.

property abort_on_unused_inputs

Set to 1 to cause the simulation to fail after its completion if there were unused parameters. (default: 0 for false) It is mainly intended for continuous integration and automated testing to check that all tests and inputs are adapted to API changes.

property always_warn_immediately

If set to 1, ImpactX immediately prints every warning message as soon as it is generated. (default: 0 for false) It is mainly intended for debug purposes, in case a simulation crashes before a global warning report can be printed.

property verbose

Controls how much information is printed to the terminal, when running ImpactX. 0 for silent, higher is more verbose. Default is 1.

property tiny_profiler

This parameter can be used to disable tiny profiling including CArena memory profiling at runtime. Default is True.

property memory_profiler

This parameter can be used to disable tiny profiler’s memory arena profiling at runtime. If `tiny_profiler is False, this parameter has no effects. Default is True.

property tiny_profiler_file

If this parameter is empty (default), the output of tiny profiling is dumped on the default out stream of AMReX. If it’s not empty, it specifies the file name for the output. Note that "/dev/null" is a special name that mean no output.

evolve()

Run the main simulation loop (deprecated, use track_particles)

track_particles()

Run the particle tracking simulation loop.

track_envelope()

Run the envelope tracking simulation loop.

Note

Our current envelope tracking implements ideal transfer maps, assuming always zero misalignments (translation or rotations). Support for misalignments and feed-down effects in envelope tracking is in development. Until then, misalignment options set on elements are silently ignored.

track_reference(ref)

Run the reference orbit tracking simulation loop.

Parameters:

ref – the reference particle (object from impactx.RefPart)

property hook

User-defined function hooks that are called, e.g, during tracking. Supported hook locations names are:

  • "before_period": before each period (e.g., turn or channel period)

  • "after_period": after each period (e.g., turn or channel period)

  • "before_element": before each element is entered

  • "after_element": after each element is exited

  • "before_slice": before each element slice

Example: Function hook that can be called before each turn (sim):

def hook_before_period(sim):
    beam = sim.beam
    turn = sim.tracking_period
    # Example: you could now manipulate elements in sim.lattice
    #          for the next turn.

sim.hook["before_period"] = hook_before_period
property tracking_step

For tracking hooks/callbacks, a global step of the simulation.

A state of internal simulation steps, increments also for space charge slice steps in elements. We start in “step 0” (initial state).

property tracking_period

For tracking hooks/callbacks, the period in the lattice (e.g., turn or channel period).

property tracking_element

For tracking hooks/callbacks, the current lattice element.

resize_mesh()

Resize the mesh domain based on the dynamic_size and related parameters.

class impactx.Config

Configuration information on ImpactX that were set at compile-time.

property have_mpi

Indicates multi-process/multi-node support via the message-passing interface (MPI). Possible values: True/False

Note

Particle beam particles are not yet dynamically load balanced. Please see the progress in issue 198.

property have_gpu

Indicates GPU support. Possible values: True/False

property gpu_backend

Indicates the available GPU support. Possible values: None, "CUDA" (for Nvidia GPUs), "HIP" (for AMD GPUs) or "SYCL" (for Intel GPUs).

property have_omp

Indicates multi-threaded CPU support via OpenMP. Possible values: True/False`

Set the environment variable OMP_NUM_THREADS to control the number of threads.

Warning

By default, OpenMP spawns as many threads as there are available physical CPU cores on a host. When MPI and OpenMP support are used at the same time, it can easily happen that one over-subscribes the available physical CPU cores. This will lead to a severe slow-down of the simulation.

By setting appropriate environment variables for OpenMP, ensure that the number of MPI processes (ranks) per node multiplied with the number of OpenMP threads is equal to the number of physical (or virtual) CPU cores. Please see our examples in the high-performance computing (HPC) on how to run efficiently in parallel environments such as supercomputers.

Particles

class impactx.ParticleContainer

Beam Particles in ImpactX.

This class stores particles, distributed over MPI ranks.

clear(keep_mass=False, keep_charge=False)

Empty the container and reset the reference particle.

Parameters:
  • keep_mass – do not reset the reference particle mass

  • keep_charge – do not reset the reference particle charge

add_n_particles(x, y, t, px, py, pt, qm, bunch_charge=None, w=None)

Add new particles to the container for fixed s.

Either the total charge (bunch_charge) or the weight of each particle (w) must be provided.

Note: This can only be used after the initialization (grids) have

been created, meaning after the call to ImpactX.init_grids() has been made in the ImpactX class.

Parameters:
  • x – positions in x

  • y – positions in y

  • t – positions as time-of-flight in c*t

  • px – momentum in x

  • py – momentum in y

  • pt – momentum in t

  • qm – charge over mass in 1/eV

  • bunch_charge – total charge within a bunch in C

  • w – weight of each particle: the macroparticle charge in units of the elementary charge e (i.e., how many real particles to represent)

ref_particle()

Access the reference particle (impactx.RefPart). Deprecated: use ref.

Returns:

return a data reference to the reference particle

Return type:

impactx.RefPart

property ref

Access the reference particle (impactx.RefPart).

set_ref_particle(refpart)

Set reference particle attributes.

Parameters:

refpart (impactx.RefPart) – a reference particle to copy all attributes from

set_bucket_length(bucket_length)

Set length of the longitudinal particle domain (e.g., length of the RF bucket in z), optionally provided for the application of particle boundary conditions.

Parameters:

bucket_length – length of the longitudinal particle domain in m

beam_moments()

Calculate beam moments at current s like the position and momentum moments of the particle distribution, as well as emittance and Twiss parameters.

Returns:

beam properties with string keywords

Return type:

dict

property store_beam_moments

In situ calculate and store the beam moments for every simulation step (default: False).

beam_moments_history()

Return the history of the beam moments on every step.

Returns:

Pandas Dataframe of beam properties, including the global reference position s

Return type:

Pandas Dataframe

reset_beam_moments_history()

Reset the history of the beam moments

min_and_max_positions()

Compute the min and max of the particle position in each dimension.

Returns:

x_min, y_min, z_min, x_max, y_max, z_max

Return type:

Tuple[float, float, float, float, float, float]

mean_and_std_positions()

Compute the mean and std of the particle position in each dimension.

Returns:

x_mean, x_std, y_mean, y_std, z_mean, z_std

Return type:

Tuple[float, float, float, float, float, float]

redistribute()

Redistribute particles in the current mesh in x, y, z.

class impactx.RefPart

This struct stores the reference particle attributes stored in impactx.ParticleContainer.

property s

integrated orbit path length, in meters

property x

horizontal position x, in meters

property y

vertical position y, in meters

property z

longitudinal position y, in meters

property t

clock time * c in meters

property px

momentum in x, normalized to mass*c, \(p_x = \gamma \beta_x\)

property py

momentum in y, normalized to mass*c, \(p_x = \gamma \beta_x\)

property pz

momentum in z, normalized to mass*c, \(p_x = \gamma \beta_x\)

property pt

energy, normalized by rest energy, \(p_t = -\gamma\)

property gamma

Read-only: Get reference particle relativistic gamma, \(\gamma = 1/\sqrt{1-\beta^2}\)

property beta

Read-only: Get reference particle relativistic beta, \(\beta = v/c\)

property beta_gamma

Read-only: Get reference particle \(\beta \cdot \gamma\)

property qm_ratio_SI

Read-only: Get reference particle charge to mass ratio (C/kg)

reset(keep_mass=False, keep_charge=False)

Reset the reference particle.

Parameters:
  • keep_mass – do not reset the reference particle mass

  • keep_charge – do not reset the reference particle charge

set_species(species_name)

Set reference particle species by name. Sets charge, mass, and gyromagnetic anomaly for a known particle species. Returns the reference particle for chaining.

Known species: electron, positron, proton, Hminus. For other species, set charge, mass, and gyromagnetic anomaly individually via set_charge_qe(), set_mass_MeV(), and set_gyromagnetic_anomaly().

Species Constants
if (species_name == "electron") {
    qe = -1.0;
    massE = m_e / MeV_inv_c2;
    g_anomaly = 0.00115965218062;
} else if (species_name == "positron") {
    qe = 1.0;
    massE = m_e / MeV_inv_c2;
    g_anomaly = 0.00115965218062;
} else if (species_name == "proton") {
    qe = 1.0;
    massE = m_p / MeV_inv_c2;
    g_anomaly = 1.7928473446;
} else if (species_name == "Hminus") {
    qe = -1.0;
    massE = (m_p + 2.0 * m_e) / MeV_inv_c2;  // neglect ~14 eV/c^2 binding energy
    g_anomaly = 1.7928473446;
}
Parameters:

species_name (str) – particle species name

Example usage:

ref = sim.beam.ref
ref.set_species("electron").set_kin_energy_MeV(2.0e3)
set_charge_qe(charge_qe)

Write-only: Set reference particle charge in (positive) elementary charges.

set_mass_MeV(massE)

Write-only: Set reference particle rest mass (MeV/c^2).

set_kin_energy_MeV(kin_energy_MeV)

Write-only: Set reference particle kinetic energy (MeV)

load_file(madx_file)

Load reference particle information from a MAD-X file.

Warning

Our MAD-X parser is under active development and provided as a preview. Please check any loaded MAD-X beams very carefully. Please report your experience and bugs on our issue tracker.

Parameters:

madx_file – file name to MAD-X file with a BEAM entry

Initial Beam Phase Space Distributions

This module provides particle beam distributions that can be used to initialize particle beams in an impactx.ParticleContainer.

Note

For additional information, consult the documentation on Beam Distribution Input. For all except the thermal distribution we allow input in two forms:

  1. Phase space ellipse axis intersections (ImpactX native)

  2. Courant-Snyder (Twiss) parameters

For the input from Twiss parameters in Python, please use the helper function twiss:

For computing Fourier coefficients from on-axis field data (used by RFCavity, SoftQuadrupole, and SoftSolenoid):

class impactx.distribution.Gaussian(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0, cutX=0.0, cutY=0.0, cutT=0.0)

A 6D Gaussian distribution, optionally with truncation. The user may specify an independent cutoff in each phase plane (x,px), (y,py), and (t,pt). The cut is performed in normalized Courant-Snyder variables corresponding to the user-supplied second moments or Twiss functions. As a result, this is equivalent to a cut corresponding to the (linearized) action in each plane. A cutoff value of 0 means no truncation (default).

Parameters:
  • lambdaX – phase space position axis intercept; for zero correlation, these are the related RMS sizes (in meters)

  • lambdaY – see lambdaX

  • lambdaT – see lambdaX

  • lambdaPx – phase space momentum axis intercept; for zero correlation, these are the related normalized RMS momenta (in radians)

  • lambdaPy – see lambdaPx

  • lambdaPt – see lambdaPx

  • muxpx – correlation length-momentum

  • muypy – see muxpx

  • mutpt – see muxpx

  • meanX – mean value of x-coordinate

  • meanY – see meanX

  • meanT – see meanX

  • meanPx – mean value of x-momentum

  • meanPy – see meanPx

  • meanPt – see meanPx

  • dispX – beam horizontal dispersion (in meters)

  • dispPx – beam horizontal dispersion derivative (dimensionless)

  • dispY – see dispX

  • dispPy – see dispPx

  • cutX – number of sigma at which to cut the distribution in (x,px) (dimensionless); 0 means no cut

  • cutY – number of sigma at which to cut the distribution in (y,py) (dimensionless); 0 means no cut

  • cutT – number of sigma at which to cut the distribution in (t,pt) (dimensionless); 0 means no cut

class impactx.distribution.Kurth4D(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A 4D Kurth distribution transversely + a uniform distribution in t + a Gaussian distribution in pt.

class impactx.distribution.Kurth6D(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A 6D Kurth distribution.

R. Kurth, Quarterly of Applied Mathematics vol. 32, pp. 325-329 (1978) C. Mitchell, K. Hwang and R. D. Ryne, IPAC2021, WEPAB248 (2021)

class impactx.distribution.KVdist(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A K-V distribution transversely + a uniform distribution in t + a Gaussian distribution in pt.

class impactx.distribution.Empty

This distribution sets all values to zero.

class impactx.distribution.Semigaussian(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A 6D Semi-Gaussian distribution (uniform in position, Gaussian in momentum).

class impactx.distribution.Triangle(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A triangle distribution for laser-plasma acceleration related applications.

A ramped, triangular current profile with a Gaussian energy spread (possibly correlated). The transverse distribution is a 4D waterbag.

class impactx.distribution.Waterbag(lambdaX, lambdaY, lambdaT, lambdaPx, lambdaPy, lambdaPt, muxpx=0.0, muypy=0.0, mutpt=0.0, meanX=0.0, meanY=0.0, meanT=0.0, meanPx=0.0, meanPy=0.0, meanPt=0.0, dispX=0.0, dispPx=0.0, dispY=0.0, dispPy=0.0)

A 6D Waterbag distribution.

class impactx.distribution.Thermal(k, kT, kT_halo, normalize, normalize_halo, halo=0.0)

A 6D stationary thermal or bithermal distribution.

Initial Beam Spin Distribution

class impactx.distribution.SpinvMF(mux, muy, muz)

A von Mises-Fisher (vMF) distribution on the unit 2-sphere.

This is used for initializing particle spin. There is a natural bijective correspondence between vMF distributions and mean (polarization) vectors.

The algorithm used here is a simplification of the algorithm described in: C. Pinzon and K. Jung, “Fast Python sampler of the von Mises Fisher distribution”, in the special case of the 2-sphere. Additional references used include:

      1. Mardia and P. E. Jupp, Directional Statistics, Wiley, 1999;

    1. Kang and H-S. Oh, “Novel sampling method for the von Mises-Fisher distribution”, Stat. and Comput. 34, 106 (2024), DOI:10.1007/s11222-024-10419-3

Parameters:
  • mux – x component of the unit vector specifying the mean direction

  • muy – y component of the unit vector specifying the mean direction

  • muz – z component of the unit vector specifying the mean direction

Lattice Elements

This module provides elements and methods for the accelerator lattice.

class impactx.elements.KnownElementsList

An iterable, list-like type of elements.

clear()

Clear the list to become empty.

extend(list)

Add a list of elements to the list.

append(element)

Add a single element to the list.

load_file(filename, nslice=1)

Load and append a lattice file from MAD-X (.madx) or PALS (e.g., .pals.yaml) formats.

Warning

Our MAD-X and PALS parsers are under active development and provided as a preview. Please check any loaded lattice files very carefully. Please report your experience and bugs on our issue tracker.

Parameters:
  • filename – filename to file with beamline elements

  • nslice – number of slices used for the application of collective effects

from_pals(pals_line, nslice=1)

Load and append a lattice from a Particle Accelerator Lattice Standard (PALS) Python Line.

Parameters:
  • pals_line – PALS Python Line with beamline elements

  • nslice – number of slices used for the application of collective effects

select(kind=None, name=None)

Filter elements by type and/or name. If both are provided, OR-based logic is applied.

Returns references to original elements, allowing modification and chaining. Chained .select(...).select(...) selections are AND-filtered.

Parameters:
  • kind – Element type(s) to filter by. Can be a string (e.g., "Drift"), regex pattern (e.g., r".*Quad"), element type (e.g., elements.Drift), or list/tuple of these.

  • name – Element name(s) to filter by. Can be a string, regex pattern, or list/tuple of these.

Return type:

impactx.elements.FilteredElementsList

Examples:

# Filter by element type
drift_elements = lattice.select(kind="Drift")
quad_elements = lattice.select(kind=elements.Quad)

# Filter by regex pattern
all_quads = lattice.select(kind=r".*Quad")  # matches Quad, ChrQuad, ExactQuad

# Filter by name
specific_elements = lattice.select(name="quad1")

# Chain filters (AND logic)
drift_named_d1 = lattice.select(kind="Drift").select(name="drift1")

# Modify original elements through references
drift_elements[0].ds = 2.0  # modifies original lattice

# delete all drifts
lattice.select(kind=r".*Drift").delete()

# replace all Quads with drift equivalents
lattice.select(kind=r".*Quad").replace_with_drifts()
get_kinds()

Get all unique element types in the lattice.

Returns:

List of unique element types (sorted by name)

Return type:

list[type]

count_by_kind(kind_pattern)

Count elements of a specific kind.

Parameters:

kind_pattern – Element kind to count. Can be string (e.g., “Drift”), regex pattern (e.g., r”.*Quad”), or element type (e.g., elements.Drift)

Returns:

Number of elements of the specified kind

Return type:

int

has_kind(kind_pattern)

Check if list contains elements of a specific kind.

Parameters:

kind_pattern – Element kind to check for. Can be string (e.g., “Drift”), regex pattern (e.g., r”.*Quad”), or element type (e.g., elements.Drift)

Returns:

True if at least one element of the specified kind exists

Return type:

bool

transfer_map(ref, order='linear', fallback_identity_map=False)

Calculate the end-to-end transfer map of the elements in the list.

Currently only the linear transfer map is implemented (order="linear"); the order parameter is reserved for future higher-order extensions. In linear mode the 6x6 map is composed element by element, using each element’s analytic per-slice linear transport map.

Collective effects like space charge, Coherent/Incoherent Synchrotron Radiation (CSR/ISR), and wakefield effects are not applied here; the returned map describes the purely linear single-particle dynamics of the design lattice.

Phase-space ordering in the returned matrix is (x, px, y, py, t, pt).

Example
#!/usr/bin/env python3
#
# Copyright 2022-2024 The ImpactX Community
#
# Authors: Axel Huebl, Chad Mitchell
# License: BSD-3-Clause-LBNL
#
# -*- coding: utf-8 -*-

import numpy as np
import pytest

from impactx import Config, RefPart, elements


def test_lattice_linear_map():
    """Calculate the linear transfer map of a lattice."""

    # Create reference particle
    ref = RefPart()
    ref.set_species("electron").set_kin_energy_MeV(1.0e3)

    # Create a valid test lattice (all elements define a linear transfer map)
    lattice = elements.KnownElementsList()
    lattice.extend(
        [
            elements.Sbend(name="bend1", ds=1.0, rc=10.0, nslice=2),
            elements.Drift(name="drift1", ds=2.0, nslice=2),
            elements.Quad(name="quad1", ds=0.5, k=1.0, nslice=2),
            elements.Drift(name="drift2", ds=1.0, nslice=2),
        ]
    )

    # Expected result (matrix multiplication)
    R_expected = np.array(
        [
            [3.74670546e-01, 2.54005827e00, 0, 0, 0, -2.34864805e-01],
            [-4.76219076e-01, -5.59489407e-01, 0, 0, 0, 3.20646253e-02],
            [0, 0, 1.64872127e00, 6.59488508e00, 0, 0],
            [0, 0, 5.21095305e-01, 2.69091188e00, 0, 0],
            [9.98334297e-02, 4.99583537e-02, 0, 0, 1.00000000e00, -1.66466013e-03],
            [0, 0, 0, 0, 0, 1.00000000e00],
        ]
    )

    if Config.precision == "SINGLE":
        atol = 1.0e-7
        rtol = 5.0e-5
    else:
        atol = 0.0
        rtol = 1.0e-8

    # Calculate Linear Transfer Map
    R = lattice.transfer_map(ref)
    assert np.allclose(R.to_numpy(), R_expected, rtol=rtol, atol=atol)

    # Check unexpected/unsupported options
    with pytest.raises(RuntimeError):
        lattice.transfer_map(ref, order="invalid")

    # Create a lattice with an element that does not define a linear transfer map
    lattice.append(elements.Programmable(name="foobar"))

    # Ensure that the calculation asserts
    with pytest.raises(RuntimeError):
        lattice.transfer_map(ref)

    # Now the user explicitly assumes that undefined maps are identity maps
    R = lattice.transfer_map(ref, fallback_identity_map=True)
    assert np.allclose(R.to_numpy(), R_expected, rtol=rtol, atol=atol)
Parameters:
  • ref – reference particle at the starting s

  • order – So far, only the calculation of linear transfer maps is supported.

  • fallback_identity_map – For elements with an undefined transfer map in the lattice, assume the identity matrix.

Returns:

The end-to-end transfer map of the lattice.

Return type:

Map6x6

map_trace(ref)

Trace the cumulative 6x6 linear transport map element by element.

The reference particle is passed by value (intentional copy); the caller’s reference particle is not modified in place. This matches the convention used by transfer_map().

This per-element trace is what twiss() consumes to transport Twiss functions through the lattice.

If you only need the final cumulative map at the lattice exit, prefer transfer_map() instead of indexing the last entry of map_trace().

Parameters:

ref – A reference particle.

Returns:

A list of dictionaries, one per lattice element plus a leading entry for the starting position. Each entry contains:

  • s – integrated path length along the reference orbit, in meters;

  • name – user-supplied element name (empty string if not named);

  • type – element type string (e.g. "Drift", "Quad", "Sbend");

  • M – cumulative 6x6 linear transport map from the start of the lattice to the exit of this element (a Map6x6 instance; call .to_numpy() for a standard C-ordered NumPy array).

The first entry always has the identity map at the starting s; the last entry contains the same map as transfer_map().

Return type:

list[dict]

to_dicts()

Serialize the lattice to a list of dictionaries.

Each element is converted to a dictionary using its to_dict() method. The resulting list can be serialized to JSON, YAML, or other formats.

Note

This transforms the buggy .to_dict() keys of ExactSbend, PlaneXYRot, PRot and ThinDipole to degrees, which by accident are written in radians. See this comment in issue #1367.

Returns:

List of element dictionaries

Return type:

list[dict]

Example:

import json
from impactx import elements

lattice = elements.KnownElementsList([
    elements.Drift(ds=1.0, name="d1"),
    elements.Quad(ds=0.5, k=2.0, name="q1"),
])

# Serialize to JSON
with open("lattice.impactx.json", "w") as f:
    json.dump(lattice.to_dicts(), f, indent=2)
from_dicts(dicts)

Load and append elements from a list of dictionaries.

Each dictionary should be in the format produced by to_dict(), containing at minimum a type key identifying the element class.

Parameters:

dicts (list[dict]) – List of element dictionaries

Example:

import json
from impactx import elements

# Load from JSON
with open("lattice.impactx.json") as f:
    data = json.load(f)

lattice = elements.KnownElementsList()
lattice.from_dicts(data)
to_py()

Generate Python code that recreates this lattice.

Returns a string containing a complete Python script with imports and a get_lattice() function that returns a KnownElementsList with all elements.

Note

Like to_dicts(), this transforms the buggy .to_dict() keys of ExactSbend, PlaneXYRot, PRot and ThinDipole from radians to degrees.

Returns:

Python source code

Return type:

str

Example:

from impactx import elements

lattice = elements.KnownElementsList([
    elements.Drift(ds=1.0, name="d1"),
    elements.Quad(ds=0.5, k=2.0, name="q1"),
])

# Generate Python code
code = lattice.to_py()
print(code)

# Save to file
with open("my_lattice.py", "w") as f:
    f.write(code)

# Later, use the generated file:
# from my_lattice import get_lattice
# lattice = get_lattice()
plot_survey(ref=None, ax=None, legend=True, legend_ncols=5)

Plot over s of all elements in the KnownElementsList.

A positive element strength denotes horizontal focusing (e.g. for quadrupoles) and bending to the right (for dipoles). In general, this depends on both the sign of the field and the sign of the charge.

Either populates the matplotlib axes in ax or creates a new axes containing the plot.

Parameters:
  • ref – A reference particle, checked for the charge sign to plot focusing/defocusing strength directions properly.

  • ax – A plotting area in matplotlib (called axes there).

  • legend – Plot a legend if true.

  • legend_ncols – Number of columns for lattice element types in the legend.

class impactx.elements.FilteredElementsList

View returned by select() on a KnownElementsList or by chained .select() on a filtered view. Indexing returns the same element objects as the full lattice; assigning to fields updates the underlying list.

All mutating operations (delete, replace_each, replace_with_drifts) rebuild the lattice using cloned elements. Existing Python references to lattice elements will then point to objects that are no longer in the lattice. If you cache element references, re-fetch them from the lattice after any mutation.

If the selection is empty, delete is a no-op and replace_* return an empty FilteredElementsList.

select(kind=None, name=None)

Narrow this view with an additional AND filter. OR logic within a single call matches select().

Parameters:
Return type:

impactx.elements.FilteredElementsList

delete()

Remove all elements in the current selection from the underlying lattice. Invalidates this view and all other live selections on that lattice. Call select() on the underlying lattice again to obtain a new view.

Return type:

None

replace_each(element, *, keep_name=True, keep_ds=False)

Replace each selected element with a copy of element. Invalidates all other live selections on the same lattice; returns a new filtered view over the same indices (the returned view is valid).

Parameters:
  • element – Element to clone at each selected index (names and ds may be overridden; see below).

  • keep_name – If true (default), copy name from each replaced element when present.

  • keep_ds – If true, copy segment length ds from each replaced element; otherwise ds comes from the template (default false).

Return type:

impactx.elements.FilteredElementsList

Examples:

# Replace quadrupoles; names kept, ds and k from template
sim.lattice.select(kind="Quad").replace_each(
    elements.Quad(name="tpl", ds=0.1, k=1.5)
)

# Same but keep ds from the replaced elements (only k from template)
sim.lattice.select(kind="Quad").replace_each(
    elements.Quad(name="tpl", ds=0.1, k=1.5), keep_ds=True
)
replace_with_drifts(*, model='match', keep_alignment=True, keep_aperture=False)

Replace each selected element with a drift of the chosen physics family. Names and ds are always taken from the replaced element. Invalidates all other live selections on the same lattice; returns a new filtered view over the same indices (the returned view is valid).

Parameters:
  • model – With "match" (default), linear elements become Drift, class names starting with Chr become ChrDrift, and class names starting with Exact become ExactDrift. With "linear", "paraxial", or "exact", every selected slot uses that drift type.

  • keep_alignment – If true (default), copy dx, dy, and rotation from each replaced element; otherwise zero them.

  • keep_aperture – If true, copy aperture_x and aperture_y from each replaced element; otherwise zero them (default false).

Return type:

impactx.elements.FilteredElementsList

Examples:

# All Quads become drifts of the matching model
sim.lattice.select(kind=r".*Quad").replace_with_drifts()

# Clear alignment errors and apertures
sim.lattice.select(kind=r".*Quad").replace_with_drifts(keep_alignment=False)
class impactx.elements.CFbend(ds, rc, k, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A combined function bending magnet. This is an ideal Sbend with a normal quadrupole field component.

Parameters:
  • ds – Segment length in m.

  • rc – Radius of curvature in m.

  • k – Quadrupole strength in m^(-2) (MADX convention) = (gradient in T/m) / (rigidity in T-m) k > 0 horizontal focusing k < 0 horizontal defocusing

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.ConstF(ds, kx, ky, kt, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A linear Constant Focusing element.

Parameters:
  • ds – Segment length in m.

  • kx – Focusing strength for x in 1/m.

  • ky – Focusing strength for y in 1/m.

  • kt – Focusing strength for t in 1/m.

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property kx

focusing x strength in 1/m

property ky

focusing y strength in 1/m

property kt

focusing t strength in 1/m

class impactx.elements.DipEdge(psi, rc, g, R=1, K0=pi**2 / 6, K1=0, K2=1, K3=1 / 6, K4=0, K5=0, K6=0, model='linear', location='entry', modify_ref_part=False, dx=0, dy=0, rotation=0, name=None)

Edge focusing associated with bend entry or exit

The model here is based on:

    1. Hwang and S. Y. Lee, PRAB 18, 122401 (2015).

as represented in the explicit, symplectic form provided in:

    1. Mitchell and K. Hwang, in Proc. NAPAC2025, TUP040, Sacramento, CA (2025).

Here, g denotes the magnetic gap, which is a length scale that sets the rate of decay of the fringe field. The values K0 - K6 denote dimensionless field integrals, describing the shape of the fringe field, as defined in eqs. (28-34) of the first reference above. In particular, K2 is the well-known fringe field parameter denoted FINT in MAD-X. The default values of the field integrals K0 - K6 are those given in eq. (52), corresponding to a tanh (i.e. logistic) field profile.

When model = "linear", the linearized map is used. This model is identical to:

      1. Brown, SLAC Report No. 75 (1982)

when expanded to first order in g/rc (gap / radius of curvature).

By comparison, note that the MAD-X DIPEDGE element uses as input the half-gap HGAP = g/2, and sets the default value FINT = 0 (while the corresponding default value of K2 is set to 1).

Note that the nonlinear model includes a nonzero horizontal translation (depending on the field integral values) that is present even for a particle that begins on the ideal “hard-edge” reference trajectory. For a beam, this will result in a centroid offset that will produce centroid oscillations in the downstream beamline. In practice, this can be avoided by aligning the downstream elements with the true horizontal position (after including the effect of the fringe field). To model this correction, we allow two options in the dipedge model:

  • the option modify_ref_part = False (default), in which the shift due to the fringe field is applied to each beam particle phase space vector but not to the reference particle phase space vector –

this model makes sense if the shift due to the fringe field is not considered in the baseline design, so that downstream elements are aligned with the “idealized” reference trajectory

  • the option modify_ref_part = True, in which the shift due to the fringe field is applied to the reference particle phase space vector, but not to the beam particle phase space vector –

this model makes sense if the shift due to the fringe field is considered as part of the baseline design, so that downstream elements are aligned with the “shifted” reference trajectory

Parameters:
  • psi – Pole face angle [radians]

  • rc – Radius of curvature [m]

  • g – Gap parameter [m]

  • R – Length scale used in fringe field integrals [m]

  • K0 – Fringe field integral [unitless]

  • K1 – Fringe field integral [unitless]

  • K2 – Fringe field integral [unitless]

  • K3 – Fringe field integral [unitless]

  • K4 – Fringe field integral [unitless]

  • K5 – Fringe field integral [unitless]

  • K6 – Fringe field integral [unitless]

  • model – the fringe field model: linear (default) or nonlinear

  • location – the fringe field edge location: entry (default) or exit

  • modify_ref_part – apply fringe field to the reference particle True or False (default)

  • dx – horizontal translation error [m]

  • dy – vertical translation error [m]

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.Drift(ds, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A drift.

Parameters:
  • ds – Segment length in m

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.ChrDrift(ds, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A drift with chromatic effects included. The Hamiltonian is expanded through second order in the transverse variables (x,px,y,py), with the exact pt dependence retained.

Parameters:
  • ds – Segment length in m

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.ExactDrift(ds, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A drift using the exact nonlinear transfer map.

Parameters:
  • ds – Segment length in m

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.Kicker(xkick, ykick, unit='dimensionless', dx=0, dy=0, rotation=0, name=None)

A thin transverse kicker.

Parameters:
  • xkick – horizontal kick strength (dimensionless OR T-m)

  • ykick – vertical kick strength (dimensionless OR T-m)

  • unit – specification of units ("dimensionless" in units of the magnetic rigidity of the reference particle or "T-m")

  • name – an optional name for the element

class impactx.elements.LinearMap(R, dx=0, dy=0, rotation=0, name=None)

A custom, linear transport matrix.

The matrix elements \(R(i,j)\) are indexed beginning with 1, so that \(i,j=1,2,3,4,5,6\).

The matrix \(R\) multiplies the phase space vector \((x,px,y,py,t,pt)\), where coordinates \((x,y,t)\) have units of m and momenta \((px,py,pt)\) are dimensionless. So, for example, \(R(1,1)\) is dimensionless, and \(R(1,2)\) has units of m.

The internal tracking methods used by ImpactX are symplectic. However, if a user-defined linear map \(R\) is provided, it is up to the user to ensure that the matrix \(R\) is symplectic. Otherwise, this condition may be violated.

Parameters:
  • R – a linear transport map to multiply with the phase space vector \((x,px,y,py,t,pt)\).

  • ds – length associated with a user-defined linear element (defaults to 0), in m

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.Multipole(multipole, K_normal, K_skew, dx=0, dy=0, rotation=0, name=None)

A general thin multipole element.

Parameters:
  • multipole – index m (m=1 dipole, m=2 quadrupole, m=3 sextupole etc.)

  • K_normal – Integrated normal multipole coefficient (meter^(-m+1)) = ds * 1/(magnetic rigidity in T-m) * (derivative of order \(m-1\) of \(B_y\) with respect to \(x\))

  • K_skew – Integrated skew multipole coefficient (meter^(-m+1))

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.ExactCFbend(ds, k_normal, k_skew, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, int_order=2, mapsteps=5, nslice=1, name=None)

A thick combined-function dipole magnet using the exact relativistic Hamiltonian, including all kinematic nonlinearities. The user must provide arrays containing normal and skew multipole coefficients, which can be specified up to decapole order. The multipole coefficients are defined in the curvilinear coordinate system defined by the nominal reference trajectory. For definitions of the coordinate system and (curvilinear) multipole coefficients we follow:

  1. Zolkin, Phys. Rev. Accel. Beams 20, 043501 (2017), DOI:10.1103/PhysRevAccelBeams.20.043501

The coefficients must appear in the following sequence:

dipole, quadrupole, sextupole, octupole, etc…

Particle tracking is performed using symplectic integration based on the Hamiltonian splitting \(H = H_1 + H_2\). Here \(H_1\) is the exact nonlinear Hamiltonian for a sector bend (including the kinematic square root), and \(H_2\) is the term containing the vector potential, which is a superposition of multipole contributions.

The vector potential is obtained from Table XI of the above-cited reference.

Parameters:
  • ds – Segment length in m.

  • k_normal – Array of normal multipole coefficients (in meter^(-m) OR in T/meter^(m-1) for m=1,2,3,..)

  • k_skew – Array of skew multipole coefficients (in meter^(-m) OR in T/meter^(m-1) for m=1,2,3,…)

  • unit – specification of units for multipole coefficients (by default, these are normalized by magnetic rigidity)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • int_order – the order used for symplectic integration (2, 4, or 6)

  • mapsteps – number of integration steps per slice used for symplectic integration

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

The arrays k_normal and k_skew must have the same number of elements.

property unit

unit specification for multipole coefficients

property int_order

the order used for symplectic integration (2, 4, or 6)

property mapsteps

number of integration steps per slice used for symplectic integration

class impactx.elements.ExactMultipole(ds, k_normal, k_skew, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, int_order=2, mapsteps=5, nslice=1, name=None)

A thick Multipole magnet using the exact relativistic Hamiltonian, including all kinematic nonlinearities. The user must provide arrays containing normal and skew multipole coefficients, which can be specified up to arbitrarily high order. The fields are assumed to be uniform along the longitudinal beamline coordinate (hard-edge model). The coefficients must appear in the following sequence:

dipole, quadrupole, sextupole, octupole, etc…

(Note: Dipole coefficients are currently ignored, and will be supported in a separate combined-function dipole element.)

Particle tracking is performed using symplectic integration based on the Hamiltonian splitting H = H_1 + H_2. Here H_1 is the nonlinear Hamiltonian for a drift (including the kinematic square root), and H_2 is the term containing the vector potential, which is a superposition of multipole contributions.

Parameters:
  • ds – Segment length in m.

  • k_normal – Array of normal multipole coefficients (in meter^(-m) OR in T/meter^(m-1) for m=1,2,3,..)

  • k_skew – Array of skew multipole coefficients (in meter^(-m) OR in T/meter^(m-1) for m=1,2,3,…)

  • unit – specification of units for multipole coefficients (by default, these are normalized by magnetic rigidity)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • int_order – the order used for symplectic integration (2, 4, or 6)

  • mapsteps – number of integration steps per slice used for symplectic integration

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property unit

unit specification for multipole coefficients

property int_order

the order used for symplectic integration (2, 4, or 6)

property mapsteps

number of integration steps per slice used for symplectic integration

class impactx.elements.Empty

This element does nothing.

class impactx.elements.NonlinearLens(knll, cnll, dx=0, dy=0, rotation=0, name=None)

Single short segment of the nonlinear magnetic insert element.

A thin lens associated with a single short segment of the nonlinear magnetic insert described by V. Danilov and S. Nagaitsev, PRSTAB 13, 084002 (2010), Sect. V.A. This element appears in MAD-X as type NLLENS.

Parameters:
  • knll – integrated strength of the nonlinear lens (m)

  • cnll – distance of singularities from the origin (m)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.BeamMonitor(name, backend='default', encoding='g', period_sample_intervals=1)

A beam monitor, writing all beam particles at fixed s to openPMD files.

If the same element name is used multiple times, then an output series is created with multiple outputs.

The I/O backend for openPMD data dumps. bp4/bp5 is the ADIOS2 I/O library, h5 is the HDF5 format, and json is a simple text format. json only works with serial/single-rank jobs. By default, the first available backend in the order given above is taken.

openPMD iteration encoding determines if multiple files are created for individual output steps or not. Variable based is an experimental feature with ADIOS2.

Parameters:
  • name – name of the series

  • backend – I/O backend, e.g., bp, h5, json

  • encoding – openPMD iteration encoding: (v)ariable based, (f)ile based, (g)roup based (default)

  • period_sample_intervals – for periodic lattice, only output every Nth period (turn)

property name

name of the series

property nonlinear_lens_invariants

Compute and output the invariants H and I within the nonlinear magnetic insert element

property alpha

Twiss alpha of the bare linear lattice at the location of output for the nonlinear IOTA invariants H and I. Horizontal and vertical values must be equal.

property beta

Twiss beta (in meters) of the bare linear lattice at the location of output for the nonlinear IOTA invariants H and I. Horizontal and vertical values must be equal.

property tn

Dimensionless strength of the IOTA nonlinear magnetic insert element used for computing H and I.

property cn

Scale factor (in meters^(1/2)) of the IOTA nonlinear magnetic insert element used for computing H and I.

class impactx.elements.Source(distribution, openpmd_path, name)

A particle source. Currently, this only supports openPMD files from our impactx.elements.BeamMonitor

Parameters:
  • distribution – Distribution type of particles in the source. currently, only "openPMD" is supported

  • openpmd_path – path to the openPMD series

  • active_once – Inject particles only for the first lattice period. Default: True

  • name – an optional name for the element

class impactx.elements.Programmable(ds=0.0, nslice=1, name=None)

A programmable beam optics element.

This element can be programmed to receive callback hooks into Python functions.

Parameters:
  • ds – Segment length in m.

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property push

This is a function hook for pushing the whole particle container. Either this function is implemented or beam_particles and ref_particle are needed. This accepts a function or lambda with the following arguments:

user_defined_function(pc: impactx.ParticleContainer, step: int)

This function is called for the particle container as it passes through the element. Note that the reference particle must be updated before the beam particles are pushed.

property beam_particles

This is a function hook for pushing all beam particles. This accepts a function or lambda with the following arguments:

user_defined_beam_function(pti: impactx.ImpactXParIter, refpart: impactx.RefPart)

This function is called repeatedly for all particle tiles or boxes in the beam particle container. Particles can be pushed and are relative to the reference particle

property ref_particle

This is a function hook for pushing the reference particle. This accepts a function or lambda with the following argument:

user_defined_refpart_function(refpart: impactx.RefPart)

This function is called for the reference particle as it passes through the element. The reference particle is updated before the beam particles are pushed.

class impactx.elements.Quad(ds, k, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A Quadrupole magnet.

Parameters:
  • ds – Segment length in m.

  • k – Quadrupole strength in m^(-2) (MADX convention) = (gradient in T/m) / (rigidity in T-m) k > 0 horizontal focusing k < 0 horizontal defocusing

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.ChrQuad(ds, k, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

A Quadrupole magnet, with chromatic effects included. The Hamiltonian is expanded through second order in the transverse variables (x,px,y,py), with the exact pt dependence retained.

Parameters:
  • ds – Segment length in m.

  • k

    Quadrupole strength in m^(-2) (MADX convention, if unit = 0)

    = (gradient in T/m) / (rigidity in T-m)

    OR Quadrupole strength in T/m (MaryLie convention, if unit = 1)

    k > 0 horizontal focusing k < 0 horizontal defocusing

  • unit – specification of units for quadrupole field strength

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property k

quadrupole strength in 1/m^2 (or T/m)

property unit

unit specification for quad strength

class impactx.elements.ExactQuad(ds, k, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, int_order=2, mapsteps=5, nslice=1, name=None)

A Quadrupole magnet using the exact relativistic Hamiltonian, including all kinematic nonlinearities. Particle tracking is performed using symplectic integration based on the Hamiltonian splitting H = H_1 + H_2. Here H_1 is the Hamiltonian for a linear quadrupole (containing all terms quadratic in the phase space variables), and H_2 is the remainder (including the kinematic square root). This suggested splitting appears for example in:

    1. Bruhwiler et al, in Proc. of EPAC 98, pp. 1171-1173 (1998).

  1. Forest, J. Phys. A: Math. Gen. 39, 5321 (2006).

Parameters:
  • ds – Segment length in m.

  • k

    Quadrupole strength in m^(-2) (MADX convention, if unit = 0)

    = (gradient in T/m) / (rigidity in T-m)

    OR Quadrupole strength in T/m (MaryLie convention, if unit = 1)

    k > 0 horizontal focusing k < 0 horizontal defocusing

  • unit – specification of units for quadrupole field strength

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • int_order – the order used for symplectic integration (2, 4, or 6)

  • mapsteps – number of integration steps per slice used for symplectic integration

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property k

quadrupole strength in 1/m^2 (or T/m)

property unit

unit specification for quad strength

property int_order

the order used for symplectic integration (2, 4, or 6)

property mapsteps

number of integration steps per slice used for symplectic integration

class impactx.elements.QuadEdge(k, unit=0, flag='entry', dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, name=None)

Hard-edge nonlinear fringe field map for a Quadrupole. This is a nonlinear symplectic map (derived from a third-order Lie generator), representing the effect of quadrupole entry or exit fringe fields in the hard-edge limit. This is an explicit symplectification of the Lie map that appears in eq (28) of: E. Forest and J. Milutinovic, Nucl. Instrum. and Methods in Phys. Res. A 269, 474-482 (1988).

Parameters:
  • k

    Quadrupole strength in m^(-2) (MADX convention, if unit = 0)

    = (gradient in T/m) / (rigidity in T-m)

    OR Quadrupole strength in T/m (MaryLie convention, if unit = 1)

    k > 0 horizontal focusing k < 0 horizontal defocusing

  • unit – specification of units for quadrupole field strength

  • flag – location of edge ("entry" or "exit")

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.ChrPlasmaLens(ds, k, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

An active cylindrically symmetric plasma lens, with chromatic effects included. The Hamiltonian is expanded through second order in the transverse variables (x,px,y,py), with the exact pt dependence retained.

Parameters:
  • ds – Segment length in m.

  • k – focusing strength in m^(-2) (if unit = 0) = (azimuthal magnetic field gradient in T/m) / (rigidity in T-m) OR azimuthal magnetic field gradient in T/m (if unit = 1)

  • unit – specification of units for plasma lens focusing strength

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property k

plasma lens focusing strength in 1/m^2 (or T/m)

property unit

unit specification for plasma lens focusing strength

class impactx.elements.ChrAcc(ds, ez, bz, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

Acceleration in a uniform field Ez, with a uniform solenoidal field Bz.

The Hamiltonian is expanded through second order in the transverse variables (x,px,y,py), with the exact pt dependence retained.

Parameters:
  • ds – Segment length in m

  • ez – electric field strength in m^(-1) = (charge * electric field Ez in V/m) / (m*c^2)

  • bz – magnetic field strength in m^(-1) = (charge * magnetic field Bz in T) / (m*c)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

property ez

electric field strength in 1/m

property bz

magnetic field strength in 1/m

class impactx.elements.RFCavity(ds, escale, freq, phase, *, cos_coefficients=None, sin_coefficients=None, z=None, field_on_axis=None, ncoef=None, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, mapsteps=10, nslice=1, name=None)

A radiofrequency cavity. See Models of Soft-Edge Elements.

Provide either pre-computed Fourier coefficients (cos_coefficients, sin_coefficients) or raw on-axis field data (z, field_on_axis, ncoef), not both. When the latter is given, Fourier coefficients are computed automatically using impactx.fourier_coefficients().

The units used for the on-axis longitudinal electric field are described in the documentation of escale below. For example, if the values used to describe the on-axis electric field (as specified in cos_coefficients, sin_coefficients, or gradient_on_axis) attain a peak on-axis value of 1, then the parameter escale, which multiplies this profile, specifies the peak value of the longitudinal electric field gradient on-axis, divided by particle rest energy.

In this case, escale has units of inverse meters.

Parameters:
  • ds – Segment length in m.

  • escale – scaling factor for on-axis RF electric field in 1/m = (peak on-axis electric field Ez in MV/m) / (particle rest energy in MeV)

  • freq – RF frequency in Hz

  • phase – RF driven phase in degrees

  • cos_coefficients – array of float cosine coefficients in Fourier expansion of on-axis electric field Ez (optional); default is a 9-cell TESLA superconducting cavity model from DOI:10.1103/PhysRevSTAB.3.092001

  • sin_coefficients – array of float sine coefficients in Fourier expansion of on-axis electric field Ez (optional); default is a 9-cell TESLA superconducting cavity model from DOI:10.1103/PhysRevSTAB.3.092001

  • z – array of longitudinal positions in m, covering the element from entry (min(z)) to exit (max(z)); the range is scaled to ds (alternative to Fourier coefficients)

  • field_on_axis – array of on-axis electric field Ez values, typically normalized to a peak absolute value of 1; multiplied by escale (alternative to Fourier coefficients)

  • ncoef – number of Fourier coefficients to compute (alternative to Fourier coefficients)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • mapsteps – number of integration steps per slice used for map and reference particle push in applied fields

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.Sbend(ds, rc, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

An ideal sector bend.

Parameters:
  • ds – Segment length in m.

  • rc – Radius of curvature in m.

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.ExactSbend(ds, phi, B=0.0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

An ideal sector bend using the exact nonlinear map. The model consists of a uniform bending field B_y with a hard edge. Pole faces are normal to the entry and exit velocity of the reference particle.

References:

      1. Bruhwiler et al, in Proc. of EPAC 98, pp. 1171-1173 (1998).

    1. Forest et al, Part. Accel. 45, pp. 65-94 (1994).

Parameters:
  • ds – Segment length in m.

  • phi – Bend angle in degrees.

  • B – Magnetic field in Tesla; when B = 0 (default), the reference bending radius is defined by r0 = length / (angle in rad), corresponding to a magnetic field of B = rigidity / r0; otherwise the reference bending radius is defined by r0 = rigidity / B.

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.Buncher(V, k, dx=0, dy=0, rotation=0)

A short RF cavity element at zero crossing for bunching (MaryLie model).

Parameters:
  • V – Normalized RF voltage drop V = Emax*L/(c*Brho)

  • k – Wavenumber of RF in 1/m

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

class impactx.elements.ShortRF(V, freq, phase=-90.0, dx=0, dy=0, rotation=0, name=None)

A short RF cavity element (MAD-X model).

Parameters:
  • V – Normalized RF voltage V = maximum energy gain/(m*c^2)

  • freq – RF frequency in Hz

  • phase – RF synchronous phase in degrees (phase = 0 corresponds to maximum energy gain, phase = -90 corresponds go zero energy gain for bunching)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.SoftSolenoid(ds, bscale, *, cos_coefficients=None, sin_coefficients=None, z=None, field_on_axis=None, ncoef=None, unit=0, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, mapsteps=10, nslice=1, name=None)

A soft-edge solenoid. See Models of Soft-Edge Elements.

Provide either pre-computed Fourier coefficients (cos_coefficients, sin_coefficients) or raw on-axis field data (z, field_on_axis, ncoef), not both. When the latter is given, Fourier coefficients are computed automatically using impactx.fourier_coefficients().

The units used for the on-axis longitudinal magnetic field data are determined by the parameter unit. For example, if the values used to describe the on-axis profile (as specified in cos_coefficients, sin_coefficients, or field_on_axis) attain a peak on-axis value of 1, then the parameter bscale, which multiplies this profile, specifies the peak value of the longitudinal magnetic field gradient on-axis. If unit=0, this is normalized by the magnetic rigidity.

Parameters:
  • ds – Segment length in m.

  • bscale – Scaling factor for on-axis magnetic field Bz in inverse meters (if unit = 0) = (magnetic field Bz in T) / (rigidity in T-m) OR Magnetic field Bz in T (SI units, if unit = 1)

  • cos_coefficients – array of float cosine coefficients in Fourier expansion of on-axis magnetic field Bz (optional); default is a thin-shell model from DOI:10.1016/J.NIMA.2022.166706

  • sin_coefficients – array of float sine coefficients in Fourier expansion of on-axis magnetic field Bz (optional); default is a thin-shell model from DOI:10.1016/J.NIMA.2022.166706

  • z – array of longitudinal positions in m, covering the element from entry (min(z)) to exit (max(z)); the range is scaled to ds (alternative to Fourier coefficients)

  • field_on_axis – array of on-axis magnetic Bz field values, typically normalized to a peak absolute value of 1; multiplied by bscale (alternative to Fourier coefficients)

  • ncoef – number of Fourier coefficients to compute (alternative to Fourier coefficients)

  • unit – specification of units for scaling of the on-axis longitudinal magnetic field

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • mapsteps – number of integration steps per slice used for map and reference particle push in applied fields

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.Sol(ds, ks, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, nslice=1, name=None)

An ideal hard-edge Solenoid magnet.

Parameters:
  • ds – Segment length in m.

  • ks – Solenoid strength in m^(-1) (MADX convention) in (magnetic field Bz in T) / (rigidity in T-m)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.PRot(phi_in, phi_out, name=None)

Exact map for a pole-face rotation in the x-z plane.

Parameters:
  • phi_in – angle of the reference particle with respect to the longitudinal (z) axis in the original frame in degrees

  • phi_out – angle of the reference particle with respect to the longitudinal (z) axis in the rotated frame in degrees

  • name – an optional name for the element

class impactx.elements.PlaneXYRot(angle, dx=0, dy=0, rotation=0, name=None)

Map for a transverse rotation in the x-y plane (i.e., about the reference velocity vector).

Parameters:
  • angle – nominal angle of rotation in the x-y plane, in degrees

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

class impactx.elements.Aperture(aperture_x, aperture_y, repeat_x, repeat_y, shift_odd_x, shape='rectangular', dx=0, dy=0, rotation=0, name=None)

A thin collimator element, applying a transverse aperture boundary.

Parameters:
  • aperture_x – horizontal half-aperture (rectangular or elliptical) in m

  • aperture_y – vertical half-aperture (rectangular or elliptical) in m

  • repeat_x – horizontal period for repeated aperture masking (inactive by default) (meter)

  • repeat_y – vertical period for repeated aperture masking (inactive by default) (meter)

  • shift_odd_x – for hexagonal/triangular mask patterns: horizontal shift of every 2nd (odd) vertical period by repeat_x / 2. Use alignment offsets dx,dy to move whole mask as needed.

  • shape – aperture boundary shape: "rectangular" (default) or "elliptical"

  • action – aperture domain action: "transmit" (default) or "absorb"

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

property shape

aperture type (rectangular, elliptical)

property action

aperture type (transmit, absorb)

property xmax

maximum horizontal coordinate

property ymax

maximum vertical coordinate

class impactx.elements.PolygonAperture(vertices_x, vertices_y, min_radius2=0.0, repeat_x, repeat_y, shift_odd_x, action="transmit", dx=0, dy=0, rotation=0, name=None)

This element defines a thin collimator element applying a transverse polygon aperture boundary defined by \((x,y)\) coordinates and optional radius below which all particles are transmitted. The vertices must define a closed curve and be ordered in the counter-clockwise direction. The first and last vertices must be identical.

Parameters:
  • vertices_x – sequence of aperture boundary \(x\) coordinates in m

  • vertices_y – sequence of aperture boundary \(y\) coordinates in m

  • min_radius2 – radius-squared of a circle fully inscribed by the polygon aperture (default 0) (meters-squared)

  • repeat_x – horizontal period for repeated aperture masking (inactive by default) (meter)

  • repeat_y – vertical period for repeated aperture masking (inactive by default) (meter)

  • shift_odd_x – for hexagonal/triangular mask patterns: horizontal shift of every 2nd (odd) vertical period by repeat_x / 2. Use alignment offsets dx,dy to move whole mask as needed.

  • action – aperture domain action: "transmit" (default) or "absorb"

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

property min_radius2

radius-squared of a fully inscribed circle. Particles with radius-squared less than this value are transmitted by the aperture and the polygon calculation is skipped.

aperture type (transmit, absorb)

class impactx.elements.SoftQuadrupole(ds, gscale, *, cos_coefficients=None, sin_coefficients=None, z=None, gradient_on_axis=None, ncoef=None, dx=0, dy=0, rotation=0, aperture_x=0, aperture_y=0, mapsteps=10, nslice=1, name=None)

A soft-edge quadrupole. See Models of Soft-Edge Elements.

Provide either pre-computed Fourier coefficients (cos_coefficients, sin_coefficients) or raw on-axis field/gradient data (z, gradient_on_axis, ncoef), not both. When the latter is given, Fourier coefficients are computed automatically using impactx.fourier_coefficients().

The units used for the on-axis quadrupole gradient are the same as those used for the quadrupole strength k in the element Quad. For example, if the values used to describe the on-axis profile (as specified in cos_coefficients, sin_coefficients, or gradient_on_axis) attain a peak on-axis value of 1, then the parameter gscale, which multiplies this profile, specifies the peak value of the quadrupole field gradient on-axis, divided by the magnetic rigidity.

In this case, gscale has units of inverse meters squared.

Parameters:
  • ds – Segment length in m.

  • gscale – Scaling factor for on-axis field gradient in inverse meters squared.

  • cos_coefficients – array of float cosine coefficients in Fourier expansion of on-axis field gradient dBy/dx (optional); default is a tanh fringe field model based on http://www.physics.umd.edu/dsat/docs/MaryLieMan.pdf

  • sin_coefficients – array of float sine coefficients in Fourier expansion of on-axis field gradient dBy/dx (optional); default is a tanh fringe field model based on http://www.physics.umd.edu/dsat/docs/MaryLieMan.pdf

  • z – array of longitudinal positions in m, covering the element from entry (min(z)) to exit (max(z)); the range is scaled to ds (alternative to Fourier coefficients)

  • gradient_on_axis – array of on-axis field gradient dBy/dx values, typically normalized to a peak absolute value of 1; multiplied by gscale (alternative to Fourier coefficients)

  • ncoef – number of Fourier coefficients to compute (alternative to Fourier coefficients)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • aperture_x – horizontal half-aperture (elliptical) in m

  • aperture_y – vertical half-aperture (elliptical) in m

  • mapsteps – number of integration steps per slice used for map and reference particle push in applied fields

  • nslice – number of slices used for the application of space charge

  • name – an optional name for the element

class impactx.elements.SpinMap(v=0, A=0, dx=0, dy=0, rotation=0, name=None)

A custom, user-specified spin map that acts on the spin 3-vector \((s_x,s_y,s_z)\). Spin maps are specified in the Lie-algebraic form:

\[\vec{s}_f = M(\zeta)\vec{s}_i,\quad\quad M(\zeta)=e^{v\cdot L}e^{A\Delta\zeta\cdot L}.\]

Here \(v\) is a 3-vector that defines the axis and angle of rotation at the phase space design point, and \(A\) is a 3x6 matrix that defines the spin-orbit coupling for particles not on the design point. Also, \(\Delta\zeta=(x,p_x,y,p_y,t,p_t)\) denotes the 6-vector of phase space variables as deviations from the design orbit. The quantities \(L_x\), \(L_y\), and \(L_z\) are standard 3x3 matrices that define a basis for the Lie algebra \(so(3)\).

The vector components \(v(i)\) and the matrix elements \(A(i,j)\) are indexed beginning with 1, so that \(i=1,2,3\) and \(j=1,2,3,4,5,6\). The vector \(v\) and the matrix \(A\) are defaulted to zero, so only entries that differ from zero need to be specified.

The matrix \(A\) multiplies the phase space vector \((x,p_x,y,p_y,t,p_t)\), where coordinates \((x,y,t)\) have units of m and momenta \((p_x,p_y,p_t)\) are dimensionless. The three components output are dimensionless. So, for example, \(A(1,1)\) has units of 1/m, and \(A(1,2)\) is dimensionless. All three components of \(v\) are dimensionless.

Parameters:
  • v – a 1-indexed, 3x1, axis-angle vector that defines the spin rotation at the phase space design point

  • R – a 1-indexed, 3x6, spin-orbit coupling matrix to multiply with the phase space vector \((x,p_x,y,p_y,t,p_t)\) that defines the spin rotation for off-design particles

  • ds – length associated with a user-defined linear element (defaults to 0), in m

  • dx – horizontal translation error in m (not used, defaults to 0)

  • dy – vertical translation error in m (not used, defaults to 0)

  • rotation – rotation error in the transverse plane [degrees] (not used, defaults to 0)

  • name – an optional name for the element

class impactx.elements.ThinDipole(theta, rc, dx=0, dy=0, rotation=0, name=None)

A general thin dipole element.

Parameters:
  • theta – Bend angle (degrees)

  • rc – Effective curvature radius (meters)

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

Reference:

    1. Ripken and F. Schmidt, Thin-Lens Formalism for Tracking, CERN/SL/95-12 (AP), 1995.

class impactx.elements.TaperedPL(k, taper, unit=0, dx=0, dy=0, rotation=0, name=None)

A thin nonlinear plasma lens with transverse (horizontal) taper

\[B_x = g \left( y + \frac{xy}{D_x} \right), \quad \quad B_y = -g \left(x + \frac{x^2 + y^2}{2 D_x} \right)\]

where \(g\) is the (linear) field gradient in T/m and \(D_x\) is the targeted horizontal dispersion in m.

Parameters:
  • k

    integrated focusing strength in m^(-1) (if unit = 0)

    = (length in m) * (magnetic field gradient \(g\) in T/m) / (magnetic rigidity in T-m)

    OR integrated focusing strength in T (if unit = 1)

    = (length in m) * (magnetic field gradient \(g\) in T/m)

  • taper – horizontal taper parameter in m^(-1) = 1 / (target horizontal dispersion \(D_x\) in m)

  • unit – specification of units for plasma lens focusing strength

  • dx – horizontal translation error in m

  • dy – vertical translation error in m

  • rotation – rotation error in the transverse plane [degrees]

  • name – an optional name for the element

property k

integrated plasma lens focusing strength in 1/m (or T)

property taper

horizontal taper parameter in 1/m

property unit

unit specification for plasma lens focusing strength

Methods

Each lattice element provides a .to_dict() method, which can be used to serialize its configuration.

elements.transformation.insert_element_every_ds(list, ds, element)

Insert an element every s into an element list.

Splits up every element that is on s = N * ds for N>0.

Parameters:
  • list – element lattice list

  • ds – spacing in meters along s to add an element

  • element – the extra element to add every s

Coordinate Transformation

class impactx.TransformationDirection

Enumerated type indicating whether to transform to fixed \(s\) or fixed \(t\) coordinate system when applying impactx.coordinate_transformation.

Parameters:
  • to_fixed_t

  • to_fixed_s

impactx.coordinate_transformation(pc, direction)

Function to transform the coordinates of the particles in a particle container either to fixed \(t\) or to fixed \(s\).

Parameters:
  • pcimpactx.particle_container whose particle coordinates are to be transformed.

  • direction – enumerated type impactx.TransformationDirection, indicates whether to transform to fixed \(s\) or fixed \(t\).