Parameters: Python
This documents on how to use ImpactX as a Python script (python3 run_script.py
).
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
, or3
- 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
toFalse
.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, turnsdynamic_size
toTrue
.
- property dynamic_size
Use dynamic (
True
) resizing of the field mesh or static sizing (False
).
- property space_charge
Enable (
True
) or disable (False
) space charge calculations (default:False
).Whether to calculate space charge effects.
- property poisson_solver
The numerical solver to solve the Poisson equation when calculating space charge effects. Either
"multigrid"
(default) or"fft"
.Currently, this is a 3D solver. An additional 2D/2.5D solver will be added in the near future.
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: ignoredThe 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
) space charge 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
andCFbend
.CSR effects require the compilation flag
-DImpactX_FFT=ON
.
- property csr_bins
Enable or disable Coherent Synchrotron Radiation (CSR) calculations (default:
150
).
- 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)
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
- particle_container()
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 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 is1
.
- evolve()
Run the main simulation loop for a number of steps.
- resize_mesh()
Resize the mesh
domain
based on thedynamic_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.
- add_n_particles(x, y, t, px, py, pt, qm, bchchg)
Add new particles to the container for fixed s.
- 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
bchchg – total charge within a bunch in C
- ref_particle()
Access the reference particle (
impactx.RefPart
).- Returns
return a data reference to the reference particle
- Return type
- set_ref_particle(refpart)
Set reference particle attributes.
- Parameters
refpart (impactx.RefPart) – a reference particle to copy all attributes from
- reduced_beam_characteristics()
Compute reduced beam characteristics 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
- 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)
- 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.
- Parameters
madx_file – file name to MAD-X file with a
BEAM
entry
Initial Beam 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:
Phase space ellipse axis intersections (ImpactX native)
Courant-Snyder (Twiss) parameters
For the input from Twiss parameters in Python, please use the helper function twiss
:
- class impactx.distribution.Gaussian(lambdax, lambday, lambdat, lambdapx, lambdapy, lambdapt, muxpx=0.0, muypy=0.0, mutpt=0.0)
A 6D Gaussian distribution.
- 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
- class impactx.distribution.Kurth4D(lambdax, lambday, lambdat, lambdapx, lambdapy, lambdapt, muxpx=0.0, muypy=0.0, mutpt=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)
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)
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)
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)
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)
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.
Lattice Elements
This module provides elements 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(madx_file, nslice=1)
Load and append an accelerator lattice description from a MAD-X file.
- Parameters
madx_file – file name to MAD-X file with beamline elements
nslice – number of slices used for the application of space charge
- class impactx.elements.CFbend(ds, rc, k, dx=0, dy=0, rotation=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]
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, 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]
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, K2, dx=0, dy=0, rotation=0, name=None)
Edge focusing associated with bend entry or exit
This model assumes a first-order effect of nonzero gap. Here we use the linear fringe field map, given to first order in g/rc (gap / radius of curvature).
References:
Brown, SLAC Report No. 75 (1982).
Hwang and S. Y. Lee, PRAB 18, 122401 (2015).
- Parameters
psi – Pole face angle in rad
rc – Radius of curvature in m
g – Gap parameter in m
K2 – Fringe field integral (unitless)
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.Drift(ds, dx=0, dy=0, rotation=0, nslice=1, name=None)
A drift.
- Parameters
ds – Segment length 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, 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]
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, 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]
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.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 (1/meter^m)
K_skew – Integrated skew multipole coefficient (1/meter^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.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.
bp
is the ADIOS2 I/O library,h5
is the HDF5 format, andjson
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.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
andref_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, 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]
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, 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]
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.ChrPlasmaLens(ds, k, unit=0, dx=0, dy=0, rotation=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]
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, 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]
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, sin_coefficients, dx=0, dy=0, rotation=0, mapsteps=1, nslice=1, name=None)
A radiofrequency cavity.
- 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.092001sin_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.092001dx – horizontal translation error in m
dy – vertical translation error in m
rotation – rotation error in the transverse plane [degrees]
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, 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]
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, 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:
Bruhwiler et al, in Proc. of EPAC 98, pp. 1171-1173 (1998).
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]
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, sin_coefficients, unit=0, dx=0, dy=0, rotation=0, mapsteps=1, nslice=1, name=None)
A soft-edge solenoid.
- 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.166706sin_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.166706unit – 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]
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, 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]
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.Aperture(xmax, ymax, shape='rectangular', dx=0, dy=0, rotation=0, name=None)
A thin collimator element, applying a transverse aperture boundary.
- Parameters
xmax – maximum allowed value of the horizontal coordinate (meter)
ymax – maximum allowed value of the vertical coordinate (meter)
shape – aperture boundary shape:
"rectangular"
(default) or"elliptical"
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 xmax
maximum horizontal coordinate
- property ymax
maximum vertical coordinate
- class impactx.elements.SoftQuadrupole(ds, gscale, cos_coefficients, sin_coefficients, dx=0, dy=0, rotation=0, mapsteps=1, nslice=1, name=None)
A soft-edge quadrupole.
- Parameters
ds – Segment length in m.
gscale – Scaling factor for on-axis field gradient in inverse meters
cos_coefficients – array of
float
cosine coefficients in Fourier expansion of on-axis field gradient (optional); default is a tanh fringe field model based on http://www.physics.umd.edu/dsat/docs/MaryLieMan.pdfsin_coefficients – array of
float
sine coefficients in Fourier expansion of on-axis field gradient (optional); default is a tanh fringe field model based on http://www.physics.umd.edu/dsat/docs/MaryLieMan.pdfdx – horizontal translation error in m
dy – vertical translation error in m
rotation – rotation error in the transverse plane [degrees]
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.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:
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
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
pc –
impactx.particle_container
whose particle coordinates are to be transformed.direction – enumerated type
impactx.TransformationDirection
, indicates whether to transform to fixed \(s\) or fixed \(t\).