Add New Beamline Elements
In ImpactX, one can easily add new beamline elements as a user. There are multiple ways to add new elements to ImpactX, you can pick the one that fits your needs best.
The workflows described here apply both for thin kicks or thick elements. Thick elements can also use soft-edged fringe fields (see existing soft-edged elements for implementation details).
Linear Map
A custom linear element can be provided by specifying the 6x6 linear transport matrix \(R\) as an input. See the FODO cell example for a demonstration of the custom linear element.
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,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. So, for example, \(R(1,1)\) is dimensionless, and \(R(1,2)\) has units of m.
Note
If a user-provided linear map is used, it is up to the user to ensure that the 6x6 transport matrix is symplectic. If a more general form of user-defined transport is needed, the Python Programmable Element and the C++ Element provide a more general approach.
Python Programmable Element
Using the ImpactX Python interface, a custom element named impactx.elements.Programmable can be defined to advance particles using NumPy, CuPy, Numba, PyTorch or any other compatible Python library.
The Programmable element can implement a custom element in two ways:
Push the whole container, by assigning a
pushfunction orPush the reference particle and beam particles in two individual functions (
beam_particlesandref_particle).
Per ImpactX convention, the reference particle is updated before the beam particles are pushed.
Detailed examples that show usage of the programmable element are:
FODO cell: implements a user-defined drift
15 stage laser-plasma accelerator: implements a user-defined LPA accelerator element using a neural network surrogate via PyTorch
Detailed particle computing interfaces are presented in the pyAMReX examples.
C++ Element
Adding a new beamline element directly to the C++ code base of ImpactX is straight forward and described in the following.
We store all beamline elements under src/elements/.
Let’s take a look at an example, the Drift implementation. To simplify the logic, we use so-called mixin classes, which provide commonly used logic for parallelization, thin/thick elements, alignment error support, etc.
struct Drift
: public mixin::Named,
public mixin::BeamOptic<Drift>,
public mixin::LinearTransport<Drift>,
public mixin::Thick,
public mixin::Alignment,
public mixin::PipeAperture,
public mixin::NoFinalize,
public amrex::simd::Vectorized<amrex::simd::native_simd_size_particlereal>
{
static constexpr auto type = "Drift";
After this brief boilerplate, our beamline elements implement these parts:
a constructor: storing element options
for particle tracking:
a single-particle operator: pushing the beam particles
a reference-particle operator: pushing the reference particle
for envelope tracking: a linear transport map
As a last step, we expose our C++ beamline elements to Python in src/python/elements.cpp.
Pull requests that added a new element and can be taken as examples are:
Chromatic Elements for Drift, Quad, Uniform Focusing+Solenoid
other pull requests under the component: elements label