FODO Channel
A 300m channel of 100 stable FODO cells (3m each) with a zero-current phase advance of 67.8 degrees.
The matched Twiss parameters at entry are:
\(\beta_\mathrm{x} = 2.82161941\) m
\(\alpha_\mathrm{x} = -1.59050035\)
\(\beta_\mathrm{y} = 2.82161941\) m
\(\alpha_\mathrm{y} = 1.59050035\)
We use a 2 GeV electron beam with initial unnormalized rms emittance of 2 nm.
The second moments of the particle distribution after the FODO cell should coincide with the second moments of the particle distribution before the FODO cell, to within the level expected due to noise due to statistical sampling.
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.
This test also demonstrates the period_sample_intervals capability of our beam monitor diagnostics, only creating output every 10th FODO cell
Run
This example can be run either as:
Python script:
python3 run_fodo.pyorImpactX executable using an input file:
impactx input_fodo.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-2024 ImpactX contributors
# Authors: Axel Huebl, Chad Mitchell, Marco Garten
# License: BSD-3-Clause-LBNL
#
# -*- coding: utf-8 -*-
from impactx import ImpactX, distribution, elements, twiss
sim = ImpactX()
# set numerical parameters and IO control
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 2 nm
kin_energy_MeV = 2.0e3 # reference energy
bunch_charge_C = 1.0e-9 # used with 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(
**twiss(
beta_x=2.8216194100262637,
beta_y=2.8216194100262637,
beta_t=0.5,
emitt_x=2e-09,
emitt_y=2e-09,
emitt_t=2e-06,
alpha_x=-1.5905003499999992,
alpha_y=1.5905003499999992,
alpha_t=0.0,
)
)
sim.add_particles(bunch_charge_C, distr, npart)
# add beam diagnostics
monitor = elements.BeamMonitor("monitor", backend="h5", period_sample_intervals=10)
# design the accelerator lattice)
ns = 5 # number of slices per ds in the element
fodo = [
monitor,
elements.Drift(ds=0.25, nslice=ns),
elements.Quad(ds=1.0, k=1.0, nslice=ns),
elements.Drift(ds=0.5, nslice=ns),
elements.Quad(ds=1.0, k=-1.0, nslice=ns),
elements.Drift(ds=0.25, nslice=ns),
]
# assign a fodo segment
sim.lattice.extend(fodo)
# FODO channel of 101 periods
sim.periods = 101
# run simulation
sim.track_particles()
# 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_from_twiss
beam.alphaX = -1.5905003499999992
beam.alphaY = -beam.alphaX
beam.alphaT = 0.0
beam.betaX = 2.8216194100262637
beam.betaY = beam.betaX
beam.betaT = 0.5
beam.emittX = 2e-09
beam.emittY = beam.emittX
beam.emittT = 2e-06
###############################################################################
# Beamline: lattice elements and segments
###############################################################################
lattice.elements = monitor drift1 quad1 drift2 quad2 drift3
lattice.nslice = 5
lattice.periods = 101 # FODO channel of 101 periods
monitor.type = beam_monitor
monitor.period_sample_intervals = 10
monitor.backend = h5
drift1.type = drift
drift1.ds = 0.25
quad1.type = quad
quad1.ds = 1.0
quad1.k = 1.0
drift2.type = drift
drift2.ds = 0.5
quad2.type = quad
quad2.ds = 1.0
quad2.k = -1.0
drift3.type = drift
drift3.ds = 0.25
###############################################################################
# Algorithms
###############################################################################
algo.space_charge = false
###############################################################################
# Diagnostics
###############################################################################
diag.slice_step_diagnostics = false
Analyze
We run the following script to analyze correctness:
Script analysis_fodo.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 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
"""
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
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)
# 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)
# compare beamline length: 300m
assert np.isclose(
300.0, series.iterations[last_step].particles["beam"].get_attribute("z_ref")
)
# compare beam monitor outputs: 10 (every 10th FODO element + 1)
assert len(series.iterations) == 11
print("Initial Beam:")
sigx, sigy, sigt, emittance_x, emittance_y, emittance_t = get_moments(initial)
print(f" sigx={sigx:e} sigy={sigy:e} sigt={sigt:e}")
print(
f" emittance_x={emittance_x:e} emittance_y={emittance_y:e} emittance_t={emittance_t:e}"
)
atol = 0.0 # ignored
rtol = 2.2 * 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],
[
7.5451170454175073e-005,
7.5441588239210947e-005,
9.9775878164077539e-004,
1.9959540393751392e-009,
2.0175015289132990e-009,
2.0013820193294972e-006,
],
rtol=rtol,
atol=atol,
)
print("")
print("Final Beam:")
sigx, sigy, sigt, emittance_x, emittance_y, emittance_t = get_moments(final)
print(f" sigx={sigx:e} sigy={sigy:e} sigt={sigt:e}")
print(
f" emittance_x={emittance_x:e} emittance_y={emittance_y:e} emittance_t={emittance_t:e}"
)
atol = 0.0 # ignored
rtol = 2.2 * 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],
[
7.4790118496224206e-005,
7.5357525169680140e-005,
9.9775879288128088e-004,
1.9959539836392703e-009,
2.0175014668882125e-009,
2.0013820380883801e-006,
],
rtol=rtol,
atol=atol,
)
Visualize
You can run the following script to visualize the beam evolution over time:
Script plot_fodo.py
#!/usr/bin/env python3
#
# Copyright 2022-2023 ImpactX contributors
# Authors: Axel Huebl, Chad Mitchell
# License: BSD-3-Clause-LBNL
#
import argparse
import glob
import re
import matplotlib.pyplot as plt
import openpmd_api as io
import pandas as pd
from matplotlib.ticker import MaxNLocator
from scipy.stats import 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
"""
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
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)
def read_file(file_pattern):
for filename in glob.glob(file_pattern):
df = pd.read_csv(filename, delimiter=r"\s+")
if "step" not in df.columns:
step = int(re.findall(r"[0-9]+", filename)[0])
df["step"] = step
yield df
def read_time_series(file_pattern):
"""Read in all CSV files from each MPI rank (and potentially OpenMP
thread). Concatenate into one Pandas dataframe.
Returns
-------
pandas.DataFrame
"""
return pd.concat(
read_file(file_pattern),
axis=0,
ignore_index=True,
) # .set_index('id')
# options to run this script
parser = argparse.ArgumentParser(description="Plot the FODO benchmark.")
parser.add_argument(
"--save-png", action="store_true", help="non-interactive run: save to PNGs"
)
args = parser.parse_args()
# 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()
ref_particle = read_time_series("diags/ref_particle.*")
# scaling to units
millimeter = 1.0e3 # m->mm
mrad = 1.0e3 # ImpactX uses "static units": momenta are normalized by the magnitude of the momentum of the reference particle p0: px/p0 (rad)
# mm_mrad = 1.e6
nm_rad = 1.0e9
# select a single particle by id
# particle_42 = beam[beam["id"] == 42]
# print(particle_42)
# steps & corresponding z
steps = list(series.iterations)
z = list(map(lambda step: ref_particle[ref_particle["step"] == step].z.values, steps))
# print(f"z={z}")
# beam transversal size & emittance over steps
moments = list(
map(
lambda step: (
step,
get_moments(series.iterations[step].particles["beam"].to_df()),
),
steps,
)
)
# print(moments)
emittance_x = list(map(lambda step_val: step_val[1][3] * nm_rad, moments))
emittance_y = list(map(lambda step_val: step_val[1][4] * nm_rad, moments))
# print(sigx, sigy)
# print beam transversal size over steps
f = plt.figure(figsize=(9, 4.8))
ax1 = f.gca()
im_emittance_x = ax1.plot(z, emittance_x, ".-", label=r"$\epsilon_x$")
im_emittance_y = ax1.plot(z, emittance_y, ".-", label=r"$\epsilon_y$")
ax1.legend()
ax1.set_xlabel(r"$z$ [m]")
ax1.set_ylabel(r"$\epsilon_{x,y}$ [nm]")
ax1.set_ylim([1.98, 2.03])
ax1.xaxis.set_major_locator(MaxNLocator(integer=True))
plt.tight_layout()
if args.save_png:
plt.savefig("fodo_lambda.png")
else:
plt.show()
Fig. 11 FODO transverse emittance evolution (preserved)