Test of a Transverse Kicker
This test applies two transverse momentum kicks, first in the horizontal direction (2 mrad) and then in the vertical direction (3 mrad).
We use a 2 GeV electron beam.
The second beam moments should be unchanged, but the first beam moments corresponding to \(p_x\) and \(p_y\) should change according to the size of the kick.
In this test, the initial and final values of \(\lambda_x\), \(\lambda_y\), \(\lambda_t\), \(\epsilon_x\), \(\epsilon_y\), and \(\epsilon_t\) must agree with nominal values.
Run
This example can be run either as:
Python script:
python3 run_kicker.py
orImpactX executable using an input file:
impactx input_kicker.in
For MPI-parallel runs, prefix these lines with mpiexec -n 4 ...
or srun -n 4 ...
, depending on the system.
#!/usr/bin/env python3
#
# Copyright 2022-2023 ImpactX contributors
# Authors: Ryan Sandberg, Axel Huebl, Chad Mitchell
# License: BSD-3-Clause-LBNL
#
# -*- coding: utf-8 -*-
from impactx import ImpactX, distribution, elements
sim = ImpactX()
# set numerical parameters and IO control
sim.particle_shape = 2 # B-spline order
sim.space_charge = False
# sim.diagnostics = False # benchmarking
sim.slice_step_diagnostics = True
# domain decomposition & space charge mesh
sim.init_grids()
# load a 2 GeV electron beam with an initial
# unnormalized rms emittance of nm
kin_energy_MeV = 2.0e3 # reference energy
bunch_charge_C = 1.0e-9 # used without space charge
npart = 10000 # number of macro particles
# reference particle
ref = sim.particle_container().ref_particle()
ref.set_charge_qe(-1.0).set_mass_MeV(0.510998950).set_kin_energy_MeV(kin_energy_MeV)
# particle bunch
distr = distribution.Waterbag(
lambdaX=4.0e-3,
lambdaY=4.0e-3,
lambdaT=1.0e-3,
lambdaPx=3.0e-4,
lambdaPy=3.0e-4,
lambdaPt=2.0e-3,
)
sim.add_particles(bunch_charge_C, distr, npart)
# add beam diagnostics
monitor = elements.BeamMonitor("monitor", backend="h5")
# design the accelerator lattice
kicklattice = [
monitor,
elements.Kicker(xkick=2.0e-3, ykick=0.0, units="dimensionless"),
elements.Kicker(xkick=0.0, ykick=3.0e-3, units="dimensionless"),
monitor,
]
# assign a lattice
sim.lattice.extend(kicklattice)
# run simulation
sim.evolve()
# clean shutdown
sim.finalize()
###############################################################################
# Particle Beam(s)
###############################################################################
beam.npart = 10000
beam.units = static
beam.kin_energy = 2.0e3
beam.charge = 1.0e-9
beam.particle = electron
beam.distribution = waterbag
beam.lambdaX = 4.0e-3
beam.lambdaY = 4.0e-3
beam.lambdaT = 1.0e-3
beam.lambdaPx = 3.0e-4
beam.lambdaPy = 3.0e-4
beam.lambdaPt = 2.0e-3
beam.muxpx = 0.0
beam.muypy = 0.0
beam.mutpt = 0.0
###############################################################################
# Beamline: lattice elements and segments
###############################################################################
lattice.elements = monitor hkick vkick monitor
monitor.type = beam_monitor
monitor.backend = h5
hkick.type = kicker
hkick.xkick = 2.0e-3 # 2 mrad horizontal kick
hkick.ykick = 0.0
vkick.type = kicker
vkick.xkick = 0.0
vkick.ykick = 3.0e-3 # 3 mrad vertical kick
###############################################################################
# Algorithms
###############################################################################
algo.particle_shape = 2
algo.space_charge = false
Analyze
We run the following script to analyze correctness:
Script analysis_kicker.py
#!/usr/bin/env python3
#
# Copyright 2022-2023 ImpactX contributors
# Authors: Axel Huebl, Chad Mitchell
# License: BSD-3-Clause-LBNL
#
import numpy as np
import openpmd_api as io
from scipy.stats import describe, moment
def get_moments(beam):
"""Calculate standard deviations of beam position & momenta
and emittance values
Returns
-------
sigx, sigy, sigt, emittance_x, emittance_y, emittance_t, meanpx, meanpy
"""
sigx = moment(beam["position_x"], moment=2) ** 0.5 # variance -> std dev.
sigpx = moment(beam["momentum_x"], moment=2) ** 0.5
sigy = moment(beam["position_y"], moment=2) ** 0.5
sigpy = moment(beam["momentum_y"], moment=2) ** 0.5
sigt = moment(beam["position_t"], moment=2) ** 0.5
sigpt = moment(beam["momentum_t"], moment=2) ** 0.5
meanpx = describe(beam["momentum_x"]).mean
meanpy = describe(beam["momentum_y"]).mean
epstrms = beam.cov(ddof=0)
emittance_x = (sigx**2 * sigpx**2 - epstrms["position_x"]["momentum_x"] ** 2) ** 0.5
emittance_y = (sigy**2 * sigpy**2 - epstrms["position_y"]["momentum_y"] ** 2) ** 0.5
emittance_t = (sigt**2 * sigpt**2 - epstrms["position_t"]["momentum_t"] ** 2) ** 0.5
return (sigx, sigy, sigt, emittance_x, emittance_y, emittance_t, meanpx, meanpy)
# initial/final beam
series = io.Series("diags/openPMD/monitor.h5", io.Access.read_only)
last_step = list(series.iterations)[-1]
initial = series.iterations[1].particles["beam"].to_df()
final = series.iterations[last_step].particles["beam"].to_df()
# compare number of particles
num_particles = 10000
assert num_particles == len(initial)
assert num_particles == len(final)
print("Initial Beam:")
sigx, sigy, sigt, emittance_x, emittance_y, emittance_t, meanpx, meanpy = get_moments(
initial
)
print(
f" sigx={sigx:e} sigy={sigy:e} sigt={sigt:e} meanpx={meanpx:e} meanpy={meanpy:e}"
)
print(
f" emittance_x={emittance_x:e} emittance_y={emittance_y:e} emittance_t={emittance_t:e}"
)
atol = 0.0 # ignored
rtol = 5.0 * num_particles**-0.5 # from random sampling of a smooth distribution
print(f" rtol={rtol} (ignored: atol~={atol})")
assert np.allclose(
[sigx, sigy, sigt, emittance_x, emittance_y, emittance_t],
[
4.017554e-03,
4.017044e-03,
9.977588e-04,
1.197572e-06,
1.210501e-06,
2.001382e-06,
],
rtol=rtol,
atol=atol,
)
atol = rtol * emittance_x / sigx # relative to rms beam size
assert np.allclose(
[meanpx, meanpy],
[
0.0,
0.0,
],
atol=atol,
)
print("")
print("Final Beam:")
sigx, sigy, sigt, emittance_x, emittance_y, emittance_t, meanpx, meanpy = get_moments(
final
)
print(
f" sigx={sigx:e} sigy={sigy:e} sigt={sigt:e} meanpx={meanpx:e} meanpy={meanpy:e}"
)
print(
f" emittance_x={emittance_x:e} emittance_y={emittance_y:e} emittance_t={emittance_t:e}"
)
atol = 0.0 # ignored
rtol = 5.0 * num_particles**-0.5 # from random sampling of a smooth distribution
print(f" rtol={rtol} (ignored: atol~={atol})")
assert np.allclose(
[sigx, sigy, sigt, emittance_x, emittance_y, emittance_t],
[
4.017554e-03,
4.017044e-03,
9.977588e-04,
1.197572e-06,
1.210501e-06,
2.001382e-06,
],
rtol=rtol,
atol=atol,
)
atol = rtol * emittance_x / sigx # relative to rms beam size
print(f" atol~={atol}")
assert np.allclose(
[meanpx, meanpy],
[
2.0e-3,
3.0e-3,
],
atol=atol,
)