HiPACE++

Tip
There will be a virtual HiPACE++ user workshop at the 11th of July 2023 from 4pm to 7pm CET. Please register on the indico webpage in case you are interested.
HiPACE++ is a 3D open-source portable (GPU-capable) quasi-static particle-in-cell code written in C++, available on GitHub. It is a full re-writing of the DESY-LBNL legacy code HiPACE, the Highly efficient Plasma Accelerator Emulator. Its main features are:
Multiple beams and plasma species to simulation beam-driven wakefield acceleration
A laser envelope solver to simulate laser-driven wakefield acceleration
An advanced explicit field solver for increased accuracy
Diagnostics compliant with the openPMD standard
Arbitrary profiles for the beams and plasma profiles
Readers from files for the beam and laser profiles
Adaptive time step and sub-cycling
Additional physics (field ionization, binary collisions, temperature effects, radiation reaction)
HiPACE++ relies on the AMReX library, which provides for particle and field data structures.
Acknowledge HiPACE++
Please acknowledge the role that HiPACE++ played in your research by citing the main HiPACE++ reference in your publication:
S. Diederichs, C. Benedetti, A. Huebl, R. Lehe, A. Myers, A. Sinn, J.-L. Vay, W. Zhang, M. Thévenet, HiPACE++: A portable, 3D quasi-static particle-in-cell code. Computer Physics Communications. ISSN 0010-4655, Volume 278, 108421, 2022 DOI:10.1016/j.cpc.2022.108421
Installation
Build/install HiPACE++
Developers
If you are new to CMake, this short tutorial from the HEP Software foundation is the perfect place to get started with it. If you just want to use CMake to build the project, jump into sections 1. Introduction, 2. Building with CMake and 9. Finding Packages.
Dependencies
HiPACE++ depends on the following popular third party software. Please see installation instructions below in the Developers section.
a mature C++17 compiler: e.g. GCC 7, Clang 7, NVCC 11.0, MSVC 19.15 or newer
AMReX development: we automatically download and compile a copy of AMReX
openPMD-api 0.15.1+: we automatically download and compile a copy of openPMD-api
Platform-dependent, at least one of the following:
CUDA Toolkit 11.0+: for NVIDIA GPU support (see matching host-compilers)
ROCm 5.2+: for AMD GPU support
FFTW3: for CPUs (only used serially, but multi-threading supported; not needed for GPUs)
Optional dependencies include:
MPI 3.0+: for multi-node and/or multi-GPU execution
OpenMP 3.1+: for threaded CPU execution
CCache: to speed up rebuilds (needs 3.7.9+ for CUDA)
Please choose one of the installation methods below to get started:
HPC platforms
If you want to use HiPACE++ on a specific high-performance computing (HPC) systems, jump directly to our HPC system-specific documentation.
Using the Spack package manager
The dependencies can be installed via the package manager Spack (macOS/Linux):
spack env create hipace-dev
spack env activate hipace-dev
spack add adios2 # for .bp file support
spack add ccache
spack add cmake
spack add fftw
spack add hdf5 # for .h5 file support
spack add mpi
spack add pkgconfig # for fftw
# optional:
# spack add cuda
spack install
(in new terminals, re-activate the environment with spack env activate hipace-dev
again)
Note
On Linux distributions, the InstallError "OpenMPI requires both C and Fortran compilers"
can occur because the Fortran compilers are sometimes not set automatically in Spack.
To fix this, the Fortran compilers must be set manually using spack config edit compilers
(more information can be found here).
For GCC, the flags f77 : null
and fc : null
must be set to f77 : gfortran
and fc : gfortran
.
On macOS, a Fortran compiler like gfortran might be missing and must be installed by hand to fix this issue.
Using the Brew package manager
The dependencies can be installed via the package manager Homebrew (macOS/Linux):
brew update
brew install adios2 # for .bp file support
brew install ccache
brew install cmake
brew install fftw
brew install hdf5-mpi # for .h5 file support
brew install libomp
brew install pkg-config # for fftw
brew install open-mpi
Now, cmake --version
should be at version 3.15.0 or newer.
Configure your compiler
For example, using a GCC on macOS:
export CC=$(which gcc)
export CXX=$(which g++)
If you also want to select a CUDA compiler:
export CUDACXX=$(which nvcc)
export CUDAHOSTCXX=$(which g++)
Build & Test
If you have not downloaded HiPACE++ yet, please clone it from GitHub via
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # or choose your preferred path
From the base of the HiPACE++ source directory, execute:
# find dependencies & configure
cmake -S . -B build
# build using up to four threads
cmake --build build -j 4
# run tests
(cd build; ctest --output-on-failure)
Note: the from_file tests require the openPMD-api with python bindings. See
documentation of the openPMD-api for more information.
An executable HiPACE++ binary with the current compile-time options encoded in its file name will be created in bin/
.
Additionally, a symbolic link named hipace
can be found in that directory, which points to the last built HiPACE++ executable. You can inspect and modify build options after running cmake .. with either
ccmake build
or by providing arguments to the CMake call
cmake -S . -B build -D<OPTION_A>=<VALUE_A> -D<OPTION_B>=<VALUE_B>
CMake Option |
Default & Values |
Description |
|
RelWithDebInfo/Release/Debug |
Type of build, symbols & optimizations |
|
NOACC/CUDA/SYCL/HIP/OMP |
On-node, accelerated computing backend |
|
ON/OFF |
Multi-node support (message-passing) |
|
SINGLE/DOUBLE |
Floating point precision (single/double) |
|
ON/OFF |
openPMD I/O (HDF5, ADIOS2) |
|
LEAPFROG/AB5 |
Use leapfrog or fifth-order Adams-Bashforth plasma pusher |
HiPACE++ can be configured in further detail with options from AMReX, which are documented in the AMReX manual.
Developers might be interested in additional options that control dependencies of HiPACE++. By default, the most important dependencies of HiPACE++ are automatically downloaded for convenience:
CMake Option |
Default & Values |
Description |
|
None |
Path to AMReX source directory (preferred if set) |
|
|
Repository URI to pull and build AMReX from |
|
|
Repository branch for |
|
ON/OFF |
Needs a pre-installed AMReX library if set to |
|
ON/OFF (default is set to value of |
Build openPMD with MPI support, although I/O is always serial |
|
None |
Path to openPMD-api source directory (preferred if set) |
|
|
Repository URI to pull and build openPMD-api from |
|
|
Repository branch for |
|
ON/OFF |
Needs a pre-installed openPMD-api library if set to |
|
ON/OFF |
Compile AMReX multigrid solver. |
For example, one can also build against a local AMReX copy.
Assuming AMReX’ source is located in $HOME/src/amrex
, add the cmake
argument -DHiPACE_amrex_src=$HOME/src/amrex
.
Relative paths are also supported, e.g. -DHiPACE_amrex_src=../amrex
.
Or build against an AMReX feature branch of a colleague.
Assuming your colleague pushed AMReX to https://github.com/WeiqunZhang/amrex/
in a branch new-feature
then pass to cmake
the arguments: -DHiPACE_amrex_repo=https://github.com/WeiqunZhang/amrex.git -DHiPACE_amrex_branch=new-feature
.
You can speed up the install further if you pre-install these dependencies, e.g. with a package manager.
Set -DHiPACE_<dependency-name>_internal=OFF
and add installation prefix of the dependency to the environment variable CMAKE_PREFIX_PATH.
Please see the short CMake tutorial that we linked in the Developers section if this sounds new to you.
Documentation
The documentation is written at the RST format, to compile the documentation locally use
cd docs
# optional: --user
python3 -m pip install -r requirements.txt # only the first time
make html
open build/html/index.html
The last line would work on MacOS. On another platform, open the html file with your favorite browser.
HPC platforms
Specific installation instructions for a set of supercomputers can be found below. Follow the guide here instead of the generic installation routines for optimal stability and best performance.
hipace.profile
Use a hipace.profile
file to set up your software environment without colliding with other software.
Ideally, store that file directly in your $HOME/
and source it after connecting to the machine:
source $HOME/hipace.profile
We list example hipace.profile
files below, which can be used to set up HiPACE++ on various HPC systems.
HPC Machines
This section documents quick-start guides for a selection of supercomputers that HiPACE++ users are active on.
Juwels Booster @ JSC
This page only provides HiPACE++ specific instructions. For more information please visit the JSC documentation.
Log in with <yourid>@juwels-booster.fz-juelich.de
.
Running on GPU
Create a file profile.hipace
and source
it whenever you log in and want to work with HiPACE++:
# please set your project account
export proj=<your project id>
# required dependencies
module load CMake
module load GCC
module load OpenMPI
module load CUDA
module load HDF5
module load ccache # optional, accelerates recompilation
export GPUS_PER_SOCKET=2
export GPUS_PER_NODE=4
# optimize CUDA compilation for A100
export AMREX_CUDA_ARCH=8.0 # 8.0 for A100, 7.0 for V100
Install HiPACE++ (the first time, and whenever you want the latest version):
source profile.hipace
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # only the first time
cd $HOME/src/hipace
rm -rf build
cmake -S . -B build -DHiPACE_COMPUTE=CUDA
cmake --build build -j 16
You can get familiar with the HiPACE++ input file format in our Get started section, to prepare an input file that suits your needs.
You can then create your directory in your $SCRATCH_<project id>
, where you can put your input file and adapt the following submission script:
#!/bin/bash -l
#SBATCH -A $proj
#SBATCH --partition=booster
#SBATCH --nodes=2
#SBATCH --ntasks=8
#SBATCH --ntasks-per-node=4
#SBATCH --gres=gpu:4
#SBATCH --time=00:05:00
#SBATCH --job-name=hipace
#SBATCH --output=hipace-%j-%N.txt
#SBATCH --error=hipace-%j-%N.err
export OMP_NUM_THREADS=1
module load GCC
module load OpenMPI
module load CUDA
module load HDF5
srun -n 8 --cpu_bind=sockets $HOME/src/hipace/build/bin/hipace.MPI.CUDA.DP.LF inputs
and use it to submit a simulation.
Tip
Parallel simulations can be largely accelerated by using GPU-aware MPI.
To utilize GPU-aware MPI, the input parameter hipace.comms_buffer_on_gpu = 1
must be set.
Note that using GPU-aware MPI may require more GPU memory.
Running on CPU
Warning
The Juwels Booster is a GPU-accelerated supercomputer, and running on CPUs only is strongly discouraged. This section only illustrates how to efficiently run on CPU with OpenMP threading, which was tested on the Juwels Booster for practical reasons, but should apply to other supercomputers. In particular, the proposed values of OMP_PROC_BIND and OMP_PLACES give decent performance for both threaded FFTW and particle operations.
Create a file profile.hipace
and source
it whenever you log in and want to work with HiPACE++:
# please set your project account
export proj=<your project id>
# required dependencies
module load CMake
module load GCC
module load OpenMPI
module load FFTW
module load HDF5
module load ccache # optional, accelerates recompilation
Install HiPACE++ (the first time, and whenever you want the latest version):
source profile.hipace
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # only the first time
cd $HOME/src/hipace
rm -rf build
cmake -S . -B build -DHiPACE_COMPUTE=OMP
cmake --build build -j 16
You can get familiar with the HiPACE++ input file format in our Get started section, to prepare an input file that suits your needs.
You can then create your directory in your $SCRATCH_<project id>
, where you can put your input file and adapt the following submission script:
#!/bin/bash -l
#SBATCH -A $proj
#SBATCH --partition=booster
#SBATCH --nodes=1
#SBATCH --ntasks=1
#SBATCH --time=00:05:00
#SBATCH --job-name=hipace
#SBATCH --output=hipace-%j-%N.txt
#SBATCH --error=hipace-%j-%N.err
source $HOME/hipace.profile
# These options give the best performance, in particular for the threaded FFTW
export OMP_PROC_BIND=false # true false master close spread
export OMP_PLACES=cores # threads cores sockets
export OMP_NUM_THREADS=8 # Anything <= 16, depending on the problem size
srun -n 8 --cpu_bind=sockets <path/to/executable> inputs
and use it to submit a simulation.
LUMI @ CSC
This page only provides HiPACE++ specific instructions. For more information please visit the LUMI documentation.
Log in with <yourid>@lumi.csc.fi
.
Running on GPU
Create a file profile.hipace
and source
it whenever you log in and want to work with HiPACE++:
# please set your project account
export proj=<your project id>
# required dependencies
module load LUMI
module load partition/G
module load PrgEnv-amd/8.3.3
module load rocm/5.2.3
module load buildtools/22.08
module load cray-hdf5/1.12.1.5
export MPICH_GPU_SUPPORT_ENABLED=1
# optimize ROCm compilation for MI250X
export AMREX_AMD_ARCH=gfx90a
# compiler environment hints
export CC=$(which cc)
export CXX=$(which CC)
export FC=$(which ftn)
Install HiPACE++ (the first time, and whenever you want the latest version):
source profile.hipace
mkdir -p $HOME/src # only the first time, create a src directory to put the codes sources
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # only the first time
cd $HOME/src/hipace
rm -rf build
cmake -S . -B build -DHiPACE_COMPUTE=HIP -DHiPACE_openpmd_mpi=OFF
cmake --build build -j 16
You can get familiar with the HiPACE++ input file format in our Get started section, to prepare an input file that suits your needs.
You can then create your directory in /project/project_<project id>
, where you can put your input file and adapt the following submission script:
#!/bin/bash
#SBATCH -A $proj
#SBATCH -J hipace
#SBATCH -o %x-%j.out
#SBATCH -t 00:30:00
#SBATCH --partition=standard-g
#SBATCH --nodes=2
#SBATCH --ntasks-per-node=8
#SBATCH --gpus-per-node=8
export MPICH_GPU_SUPPORT_ENABLED=1
# note (12-12-22)
# this environment setting is currently needed on LUMI to work-around a
# known issue with Libfabric
#export FI_MR_CACHE_MAX_COUNT=0 # libfabric disable caching
# or, less invasive:
export FI_MR_CACHE_MONITOR=memhooks # alternative cache monitor
# note (9-2-22, OLCFDEV-1079)
# this environment setting is needed to avoid that rocFFT writes a cache in
# the home directory, which does not scale.
export ROCFFT_RTC_CACHE_PATH=/dev/null
export OMP_NUM_THREADS=1
# needed for high nz runs.
# if too many mpi messages are send, the hardware counters can overflow, see
# https://docs.nersc.gov/performance/network/
# export FI_CXI_RX_MATCH_MODE=hybrid
# setting correct CPU binding
# (see https://docs.lumi-supercomputer.eu/runjobs/scheduled-jobs/lumig-job/)
cat << EOF > select_gpu
#!/bin/bash
export ROCR_VISIBLE_DEVICES=\$SLURM_LOCALID
exec \$*
EOF
chmod +x ./select_gpu
CPU_BIND="map_cpu:49,57,17,25,1,9,33,41"
srun --cpu-bind=${CPU_BIND} ./select_gpu $HOME/src/hipace/build/bin/hipace inputs
rm -rf ./select_gpu
and use it to submit a simulation.
Tip
Parallel simulations can be largely accelerated by using GPU-aware MPI.
To utilize GPU-aware MPI, the input parameter hipace.comms_buffer_on_gpu = 1
must be set and the following flag must be passed in the job script:
export FI_MR_CACHE_MAX_COUNT=0
Note that using GPU-aware MPI may require more GPU memory.
Maxwell cluster @ DESY
This page only provides HiPACE++ specific instructions. For more information please visit the Maxwell documentation.
Create a file profile.hipace
and source
it whenever you log in and want to work with
HiPACE++:
module purge
module load maxwell gcc/9.3 openmpi/4
module load maxwell cuda/11.8
module load hdf5/1.10.6
# pick correct GPU setting (this may differ for V100 nodes)
export GPUS_PER_SOCKET=2
export GPUS_PER_NODE=4
# optimize CUDA compilation for A100
export AMREX_CUDA_ARCH=8.0 # use 8.0 for A100 or 7.0 for V100
Install HiPACE++ (the first time, and whenever you want the latest version):
source profile.hipace
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # only the first time
cd $HOME/src/hipace
rm -rf build
cmake -S . -B build -DHiPACE_COMPUTE=CUDA
cmake --build build -j 16
You can get familiar with the HiPACE++ input file format in our Get started
section, to prepare an input file that suits your needs. You can then create your directory on
BEEGFS /beegfs/desy/group/<your group>
or /beegfs/desy/user/<your username>
,
where you can put your input file and adapt the following
submission script:
#! /usr/bin/env sh
#SBATCH --partition=<partition> # mpa # maxgpu # allgpu
#SBATCH --time=01:00:00
#SBATCH --nodes=1
#SBATCH --constraint=A100&GPUx4 # A100&GPUx1
#SBATCH --job-name=HiPACE
#SBATCH --output=hipace-%j-%N.out
#SBATCH --error=hipace-%j-%N.err
export OMP_NUM_THREADS=1
module load maxwell gcc/9.3 openmpi/4 cuda/11.8
mpiexec -n 4 -npernode 4 $HOME/src/hipace/build/bin/hipace inputs
The -npernode
must be set to GPUS_PER_NODE
, otherwise not all GPUs are used correctly.
There are nodes with 4 GPUs and 1 GPU (see the Maxwell documentation on compute infrastructure.
for more details and the required constraints). Please set the value accordingly.
Tip
Parallel simulations can be largely accelerated by using GPU-aware MPI.
To utilize GPU-aware MPI, the input parameter hipace.comms_buffer_on_gpu = 1
must be set and the following flag must be passed in the job script:
export UCX_MEMTYPE_CACHE=n
Note that using GPU-aware MPI may require more GPU memory.
Perlmutter @ NERSC
This page only provides HiPACE++ specific instructions. For more information please visit the Perlmutter documentation.
Log in with ssh <yourid>@perlmutter-p1.nersc.gov
.
Building for GPU
Create a file profile.hipace
and source
it whenever you log in and want to work with HiPACE++:
# please set your project account
export proj=<your project id>_g # _g for GPU accounting
# required dependencies
module load cmake/3.22.0
module load cray-hdf5-parallel/1.12.2.3
# necessary to use CUDA-Aware MPI and run a job
export CRAY_ACCEL_TARGET=nvidia80
# optimize CUDA compilation for A100
export AMREX_CUDA_ARCH=8.0
# compiler environment hints
export CC=cc
export CXX=CC
export FC=ftn
export CUDACXX=$(which nvcc)
export CUDAHOSTCXX=CC
Download HiPACE++ from GitHub (the first time, and whenever you want the latest version):
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # or any other path you prefer
Compile the code using CMake
source profile.hipace # load the correct modules
cd $HOME/src/hipace # or where HiPACE++ is installed
rm -rf build
cmake -S . -B build -DHiPACE_COMPUTE=CUDA
cmake --build build -j 16
You can get familiar with the HiPACE++ input file format in our Get started section, to prepare an input file that suits your needs.
You can then create your directory in your $PSCRATCH
, where you can put your input file and adapt the following submission script:
#!/bin/bash -l
#SBATCH -t 01:00:00
#SBATCH -N 2
#SBATCH -J HiPACE++
# note: <proj> must end on _g
#SBATCH -A <proj>_g
#SBATCH -q regular
#SBATCH -C gpu
#SBATCH -c 32
#SBATCH --exclusive
#SBATCH --gpu-bind=none
#SBATCH --gpus-per-node=4
#SBATCH -o hipace.o%j
#SBATCH -e hipace.e%j
# path to executable and input script
EXE=$HOME/src/hipace/build/bin/hipace
INPUTS=inputs
# pin to closest NIC to GPU
export MPICH_OFI_NIC_POLICY=GPU
# for GPU-aware MPI use the first line
#HIPACE_GPU_AWARE_MPI="hipace.comms_buffer_on_gpu=1"
HIPACE_GPU_AWARE_MPI=""
# CUDA visible devices are ordered inverse to local task IDs
# Reference: nvidia-smi topo -m
srun --cpu-bind=cores bash -c "
export CUDA_VISIBLE_DEVICES=\$((3-SLURM_LOCALID));
${EXE} ${INPUTS} ${HIPACE_GPU_AWARE_MPI}" \
> output.txt
and use it to submit a simulation. Note, that this example simulation runs on 8 GPUs, since -N = 2 yields 2 nodes with 4 GPUs each.
Tip
Parallel simulations can be largely accelerated by using GPU-aware MPI.
To utilize GPU-aware MPI, the input parameter hipace.comms_buffer_on_gpu = 1
must be set (see the job script above).
Note that using GPU-aware MPI may require more GPU memory.
Spock @ OLCF
This page only provides HiPACE++ specific instructions. For more information please visit the OLCF documentation.
Running on GPU
Create a file profile.hipace.spock
and source
it whenever you log in and want to work with HiPACE++:
module load cmake/3.20.2 # cmake
module load craype-accel-amd-gfx908
module load rocm # rocm/4.5.0
module load ccache
module load ninja
export AMREX_AMD_ARCH=gfx908
export CC=$(which clang)
export CXX=$(which hipcc)
export LDFLAGS="-L${CRAYLIBS_X86_64} $(CC --cray-print-opts=libs) -lmpi"
Install HiPACE++ (the first time, and whenever you want the latest version):
source profile.hipace
git clone https://github.com/Hi-PACE/hipace.git $HOME/src/hipace # only the first time
cd $HOME/src/hipace
cmake -S . -B build -DHiPACE_COMPUTE=HIP -DAMReX_AMD_ARCH=gfx908 -DMPI_CXX_COMPILER=$(which CC) -DMPI_C_COMPILER=$(which cc) -DMPI_COMPILER_FLAGS="--cray-print-opts=all"
cmake --build build -j 6
You can get familiar with the HiPACE++ input file format in our Get started section, to prepare an input file that suits your needs.
You can then create your directory in your $SCRATCH_<project id>
, where you can put your input file and adapt the following submission script:
#!/bin/bash
#SBATCH -A aph114
#SBATCH -J hipace
#SBATCH -o %x-%j.out
#SBATCH -t 00:10:00
#SBATCH -p ecp
#SBATCH -N 1
source $HOME/profile.hipace.spock
export OMP_NUM_THREADS=1
srun -n 1 -c 1 --ntasks-per-node=1 <path/to/executable>/hipace.MPI.HIP.DP.LF inputs &> output.txt
and use it to submit a simulation.
Tip
Your HPC system is not in the list? Open an issue and we can document it together!
Usage
Input parameters
Parser
In HiPACE++ all input parameters are obtained through amrex::Parser
, making it possible to
specify input parameters with expressions and not just numbers. User constants can be defined
in the input script with my_constants
.
my_constants.ne = 1.25e24
my_constants.kp_inv = "clight / sqrt(ne * q_e^2 / (epsilon0 * m_e))"
beam.radius = "kp_inv / 2"
Thereby, the following constants are predefined:
variable |
name |
Value |
q_e |
elementary charge |
1.602176634e-19 |
m_e |
electron mass |
9.1093837015e-31 |
m_p |
proton mass |
1.67262192369e-27 |
epsilon0 |
vacuum permittivity |
8.8541878128e-12 |
mu0 |
vacuum permeability |
1.25663706212e-06 |
clight |
speed of light |
299’792’458. |
hbar |
reduced Planck constant |
1.054571817e-34 |
r_e |
classical electron radius |
2.817940326204929e-15 |
For a list of supported functions see the
AMReX documentation.
Sometimes it is necessary to use double-quotes around expressions, especially when providing them
as command line parameters. Multi-line expressions are allowed if surrounded by double-quotes.
For input parameters of type string
, my_constants
can be putinside curly braces {...}
to directly paste them into the input parameter.
If what is inside the braces is not a my_constants
it will be evaluated as an expression using the parser.
General parameters
hipace.normalized_units
(bool) optional (default 0)Whether the simulation uses the normalized unit system commonly used in wakefield acceleration, see e.g. chapter 2 of this reference. Otherwise, the code assumes SI (Système International) unit system.
random_seed
(integer) optional (default 1)Passes a seed to the AMReX random number generator. This allows for reproducibility of random events such as randomly generated beams, ionization, and collisions. Note that on GPU, since the order of operations is not ensured, providing a seed does not guarantee reproducibility to the level of machine precision.
hipace.verbose
(int) optional (default 0)Level of verbosity.
hipace.verbose = 1
, prints only the time steps, which are computed.hipace.verbose = 2
additionally prints the number of iterations in the predictor-corrector loop, as well as the B-Field error at each slice.hipace.verbose = 3
also prints the number of particles that violate the quasi-static approximation and were neglected at each slice. It prints the number of ionized particles, if ionization occurred. It also adds additional information if beams are read in from file.
hipace.do_device_synchronize
(int) optional (default 0)Level of synchronization on GPU.
hipace.do_device_synchronize = 0
, synchronization happens only when necessary.hipace.do_device_synchronize = 1
, synchronizes most functions (all that are profiled viaHIPACE_PROFILE
)hipace.do_device_synchronize = 2
additionally synchronizes low-level functions (all that are profiled viaHIPACE_DETAIL_PROFILE
)
amrex.the_arena_is_managed
(bool) optional (default 0)Whether managed memory is used. Note that large simulations sometimes only fit on a GPU if managed memory is used, but generally it is recommended to not use it.
amrex.omp_threads
(system
,nosmt
or positive integer; default isnosmt
)An integer number can be set in lieu of the
OMP_NUM_THREADS
environment variable to control the number of OpenMP threads to use for theOMP
compute backend on CPUs. By default, we use thenosmt
option, which overwrites the OpenMP default of spawning one thread per logical CPU core, and instead only spawns a number of threads equal to the number of physical CPU cores on the machine. If set, the environment variableOMP_NUM_THREADS
takes precedence oversystem
andnosmt
, but not over integer numbers set in this option.
hipace.comms_buffer_on_gpu
(bool) optional (default 0)Whether the buffers that hold the beam and the 3D laser envelope should be allocated on the GPU (device memory). By default they will be allocated on the CPU (pinned memory). Setting this option to 1 is necessary to take advantage of GPU-Enabled MPI, however for this additional enviroment variables need to be set depending on the system.
hipace.comms_buffer_max_leading_slices
(int) optional (default inf)How many slices of beam particles can be received and stored in advance.
hipace.comms_buffer_max_trailing_slices
(int) optional (default inf)How many slices of beam particles can be stored before being sent. Using
comms_buffer_max_leading_slices
andcomms_buffer_max_trailing_slices
will in principle limit the amount of asynchronousness in the parallel communication and may thus reduce performance. However it may be necessary to set these parameters to avoid all slices accumulating on a single rank that would run out of memory (out of CPU or GPU memory depending onhipace.comms_buffer_on_gpu
). If there are more time steps than ranks, these parameters must be chosen such that between all ranks there is enough capacity to store every slice to avoid a deadlock, i.e.(comms_buffer_max_leading_slices + comms_buffer_max_trailing_slices) * nranks > nslices
.
hipace.do_tiling
(bool) optional (default true)Whether to use tiling, when running on CPU. Currently, this option only affects plasma operations (gather, push and deposition). The tile size can be set with
plasmas.sort_bin_size
.
hipace.depos_order_xy
(int) optional (default 2)Transverse particle shape order. Currently, 0,1,2,3 are implemented.
hipace.depos_order_z
(int) optional (default 0)Longitudinal particle shape order. Currently, only 0 is implemented.
hipace.depos_derivative_type
(int) optional (default 2)Type of derivative used in explicit deposition. 0: analytic, 1: nodal, 2: centered
hipace.outer_depos_loop
(bool) optional (default 0)If the loop over depos_order is included in the loop over particles.
hipace.do_beam_jx_jy_deposition
(bool) optional (default 1)Using the default, the beam deposits all currents
Jx
,Jy
,Jz
. Usinghipace.do_beam_jx_jy_deposition = 0
disables the transverse current deposition of the beams.
hipace.do_beam_jz_minus_rho
(bool) optional (default 0)Whether the beam contribution to \(j_z-c\rho\) is calculated and used when solving for Psi (used to caculate the transverse fields Ex-By and Ey+Bx). if 0, this term is assumed to be 0 (a good approximation for an ultra-relativistic beam in the z direction with small transverse momentum).
hipace.interpolate_neutralizing_background
(bool) optional (default 0)Whether the neutralizing background from plasmas should be interpolated from level 0 to higher MR levels instead of depositing it on all levels.
hipace.output_input
(bool) optional (default 0)Print all input parameters before running the simulation. If a parameter is present multiple times then the last occurrence will be used. Note that this will include some default AMReX parameters.
Geometry
amr.n_cell
(3 integer)Number of cells in x, y and z. With the explicit solver (default), the number of cells in the x and y directions must be either \(2^n-1\) (common values are 511, 1023, 2047, best configuration for performance) or \(2^n\) where \(n\) is an integer. Some other values might work, like \(3 \times 2^n-1\), but use at your own risk.
amr.max_level
(integer) optional (default 0)Maximum level of mesh refinement. Currently, mesh refinement is supported up to level 2. Note, that the mesh refinement algorithm is still in active development and should be used with care.
geometry.patch_lo
(3 float)Lower end of the simulation box in x, y and z.
geometry.patch_hi
(3 float)Higher end of the simulation box in x, y and z.
geometry.is_periodic
(3 bool)Whether the boundary conditions for particles in x, y and z is periodic. Note that particles in z are always removed. This setting will most likely be changed in the near future.
mr_lev1.n_cell
(2 integer)Number of cells in x and y for level 1. The number of cells in the zeta direction is calculated from
patch_lo
andpatch_hi
.
mr_lev1.patch_lo
(3 float)Lower end of the refined grid in x, y and z.
mr_lev1.patch_hi
(3 float)Upper end of the refined grid in x, y and z.
mr_lev2.n_cell
(2 integer)Number of cells in x and y for level 2. The number of cells in the zeta direction is calculated from
patch_lo
andpatch_hi
.
mr_lev2.patch_lo
(3 float)Lower end of the refined grid in x, y and z.
mr_lev2.patch_hi
(3 float)Upper end of the refined grid in x, y and z.
Time step
max_step
(integer) optional (default 0)Maximum number of time steps. 0 means that the 0th time step will be calculated, which are the fields of the initial beams.
hipace.max_time
(float) optional (default infinity)Maximum physical time of the simulation. The
dt
of the last time step may be reduced so thatt + dt = max_time
, both for the adaptive and a fixed time step.
hipace.dt
(float or string) optional (default 0.)Time step to advance the particle beam. For adaptive time step, use
"adaptive"
.
hipace.dt_max
(float) optional (default inf)Only used if
hipace.dt = adaptive
. Upper bound of the adaptive time step: if the computed adaptive time step is is larger thandt_max
, thendt_max
is used instead. Useful when the plasma profile starts with a very low density (e.g. in the presence of a realistic density ramp), to avoid unreasonably large time steps.
hipace.nt_per_betatron
(Real) optional (default 20.)Only used when using adaptive time step (see
hipace.dt
above). Number of time steps per betatron period (of the full blowout regime). The time step is given by \(\omega_{\beta}\Delta t = 2 \pi/N\) (\(N\) isnt_per_betatron
) where \(\omega_{\beta}=\omega_p/\sqrt{2\gamma}\) with \(\omega_p\) the plasma angular frequency and \(\gamma\) is an average of Lorentz factors of the slowest particles in all beams.
hipace.adaptive_predict_step
(bool) optional (default 1)Only used when using adaptive time step (see
hipace.dt
above). If true, the current Lorentz factor and accelerating field on the beams are used to predict the (adaptive)dt
of the next time steps. This prediction is used to better estimate the betatron frequency at the beginning of the next step performed by the current rank. It improves accuracy for parallel simulations (with significant deceleration and/or z-dependent plasma profile). Note: should be on by default once good defaults are determined.
hipace.adaptive_control_phase_advance
(bool) optional (default 1)Only used when using adaptive time step (see
hipace.dt
above). If true, a test on the phase advance sets the time step so it matches the phase advance expected for a uniform plasma (to a certain tolerance). This should improve the accuracy in the presence of density gradients. Note: should be on by default once good defaults are determined.
hipace.adaptive_phase_tolerance
(Real) optional (default 4.e-4)Only used when using adaptive time step (see
hipace.dt
above) andadaptive_control_phase_advance
. Tolerance for the controlled phase advance described above (lower is more accurate, but should result in more time steps).
hipace.adaptive_phase_substeps
(int) optional (default 2000)Only used when using adaptive time step (see
hipace.dt
above) andadaptive_control_phase_advance
. Number of sub-steps in the controlled phase advance described above (higher is more accurate, but should be slower).
hipace.adaptive_threshold_uz
(Real) optional (default 2.)Only used when using adaptive time step (see
hipace.dt
above). Threshold beam momentum, below which the time step is not decreased (to avoid arbitrarily small time steps).
Field solver parameters
Two different field solvers are available to calculate the transverse magnetic fields Bx and By: an explicit solver (based on analytic integration) and a predictor-corrector loop (based on an FFT solver). In the explicit solver, the longitudinal derivative of the transverse currents is calculated explicitly, which results in a shielded Poisson equation, solved with either the internal HiPACE++ multigrid solver or the AMReX multigrid solver. The default is to use the explicit solver. We strongly recommend to use the explicit solver, because we found it to be more robust, faster to converge, and easier to use.
hipace.bxby_solver
(string) optional (default explicit)Which solver to use. Possible values:
explicit
andpredictor-corrector
.
hipace.use_small_dst
(bool) optional (default 0 or 1)Whether to use a large R2C or a small C2R fft in the dst of the Poisson solver. The small dst is quicker for simulations with \(\geq 511\) transverse grid points. The default is set accordingly.
fields.extended_solve
(bool) optional (default 0)Extends the area of the FFT Poisson solver to the ghost cells. This can reduce artifacts originating from the boundary for long simulations.
fields.open_boundary
(bool) optional (default 0)Uses a Taylor approximation of the Greens function to solve the Poisson equations with open boundary conditions. It’s recommended to use this together with
fields.extended_solve = true
andgeometry.is_periodic = false false false
. Only available with the predictor-corrector solver.
fields.do_symmetrize
(bool) optional (default 0)Symmetrizes current and charge densities transversely before the field solve. Each cell at (x, y) is averaged with cells at (-x, y), (x, -y) and (-x, -y).
Explicit solver parameters
hipace.use_amrex_mlmg
(bool) optional (default 0)Whether to use the AMReX multigrid solver. Note that this requires the compile-time option
AMReX_LINEAR_SOLVERS
to be true. Generally not recommended since it is significantly slower than the default HiPACE++ multigrid solver.
hipace.MG_tolerance_rel
(float) optional (default 1e-4)Relative error tolerance of the multigrid solvers.
hipace.MG_tolerance_abs
(float) optional (default 0.)Absolute error tolerance of the multigrid solvers.
hipace.MG_verbose
(int) optional (default 0)Level of verbosity of the the multigrid solvers.
Predictor-corrector loop parameters
hipace.predcorr_B_error_tolerance
(float) optional (default 4e-2)The tolerance of the transverse B-field error. Set to a negative value to use a fixed number of iterations.
hipace.predcorr_max_iterations
(int) optional (default 30)The maximum number of iterations in the predictor-corrector loop for single slice.
hipace.predcorr_B_mixing_factor
(float) optional (default 0.05)The mixing factor between the currently calculated B-field and the B-field of the previous iteration (or initial guess, in case of the first iteration). A higher mixing factor leads to a faster convergence, but increases the chance of divergence.
Note
In general, we recommend two different settings:
First, a fixed B-field error tolerance. This ensures the same level of convergence at each grid
point. To do so, use e.g. the default settings of hipace.predcorr_B_error_tolerance = 4e-2
,
hipace.predcorr_max_iterations = 30
, hipace.predcorr_B_mixing_factor = 0.05
.
This should almost always give reasonable results.
Second, a fixed (low) number of iterations. This is usually much faster than the fixed B-field
error, but can loose significant accuracy in special physical simulation settings. For most
settings (e.g. a standard PWFA simulation the blowout regime at a reasonable resolution) it
reproduces the same results as the fixed B-field error tolerance setting. It works very well at
high longitudinal resolution.
A good setting for the fixed number of iterations is usually given by
hipace.predcorr_B_error_tolerance = -1.
, hipace.predcorr_max_iterations = 1
,
hipace.predcorr_B_mixing_factor = 0.15
. The B-field error tolerance must be negative.
Plasma parameters
The name of all plasma species must be specified with plasmas.names = ….
Then, properties can be set per plasma species with <plasma name>.<plasma property> = ...
,
or sometimes for all plasma species at the same time with plasmas.<plasma property> = ...
.
When both are specified, the per-species value is used.
plasmas.names
(string) optional (default no_plasma)The names of the plasmas, separated by a space. To run without plasma, choose the name
no_plasma
.
<plasma name> or plasmas.density(x,y,z)
(float) optional (default 0.)The plasma density as function of x, y and z. x and y coordinates are taken from the simulation box and \(z = time \cdot c\). The density gets recalculated at the beginning of every timestep. If specified as a command line parameter, quotation marks must be added:
"<plasma name>.density(x,y,z)" = "1."
.
<plasma name> or plasmas.min_density
(float) optional (default 0)Particles with a density less than or equal to the minimal density won’t be injected. Useful for parsed functions to avoid redundant plasma particles with close to 0 weight.
<plasma name>.density_table_file
(string) optional (default “”)Alternative to
<plasma name>.density(x,y,z)
. Specify the name of a text file containing multiple densities for different positions. File syntax:<position> <density function>
for every line. If a line doesn’t start with a position it is ignored (comments can be made with #). <density function> is evaluated like<plasma name>.density(x,y,z)
. The simulation position \(time \cdot c\) is rounded up to the nearest <position> in the file to get it’s <density function> which is used for that time step.
<plasma name> or plasmas.ppc
(2 integer)The number of plasma particles per cell in x and y. Since in a quasi-static code, there is only a 2D plasma slice evolving along the longitudinal coordinate, there is no need to specify a number of particles per cell in z.
<plasma name> or plasmas.radius
(float) optional (default infinity)Radius of the plasma. Set a value to run simulations in a plasma column.
<plasma name> or plasmas.hollow_core_radius
(float) optional (default 0.)Inner radius of a hollow core plasma. The hollow core radius must be smaller than the plasma radius itself.
<plasma name> or plasmas.max_qsa_weighting_factor
(float) optional (default 35.)The maximum allowed weighting factor \(\gamma /(\psi+1)\) before particles are considered as violating the quasi-static approximation and are removed from the simulation.
<plasma name>.mass
(float) optional (default 0.)The mass of plasma particle in SI units. Use
plasma_name.mass_Da
for Dalton. Can also be set with<plasma name>.element
. Must be >0.
<plasma name>.mass_Da
(float) optional (default 0.)The mass of plasma particle in Dalton. Use
<plasma name>.mass
for SI units. Can also be set with<plasma name>.element
. Must be >0.
<plasma name>.charge
(float) optional (default 0.)The charge of a plasma particle. Can also be set with
<plasma name>.element
. The charge gets multiplied by the current ionization level.
<plasma name>.element
(string) optional (default “”)The physical element of the plasma. Sets charge, mass and, if available, the specific ionization energy of each state. Options are:
electron
,positron
,H
,D
,T
,He
,Li
,Be
,B
, ….
<plasma name>.can_ionize
(bool) optional (default 0)Whether this plasma can ionize. Can also be set to 1 by specifying
<plasma name>.ionization_product
.
<plasma name>.initial_ion_level
(int) optional (default -1)The initial ionization state of the plasma. 0 for neutral gasses. If set, the plasma charge gets multiplied by this number. If the plasma species is not ionizable, the initial ionization level is set to 1.
<plasma name>.ionization_product
(string) optional (default “”)Name of the plasma species that contains the new electrons that are produced when this plasma gets ionized. Only needed if this plasma is ionizable.
<plasma name> or plasmas.neutralize_background
(bool) optional (default 1)Whether to add a neutralizing background of immobile particles of opposite charge.
plasmas.sort_bin_size
(int) optional (default 32)Tile size for plasma current deposition, when running on CPU. When tiling is activated (
hipace.do_tiling = 1
), the current deposition is done in temporary arrays of sizesort_bin_size
(+ guard cells) that are atomic-added to the main current arrays.
<plasma name>.temperature_in_ev
(float) optional (default 0)- Initializes the plasma particles with a given temperature \(k_B T\) in eV. Using a temperature, the plasma particle momentum is normally distributed with a variance of \(k_B T /(M c^2)\) in each dimension, with \(M\) the particle mass, \(k_B\) the Boltzmann constant, and \(T\) the isotropic temperature in Kelvin.Note: Using a temperature can affect the performance since the plasma particles loose their order and thus their favorable memory access pattern. The performance can be mostly recovered by reordering the plasma particles (see
<plasma name> or plasmas.reorder_period
). Furthermore, the noise of the temperature can seed the hosing instability. The amplitude of the seeding is unphysical, because the number of macro-particles is typically orders of magnitude below the number of actual plasma electrons. Since it is often unfeasible to use a sufficient amount of plasma macro-particles per cell to suppress this numerical seed, the plasma can be symmetrized to prevent the onset of the hosing instability (see<plasma name> or plasmas.do_symmetrize
).
<plasma name> or plasmas.do_symmetrize
(bool) optional (default 0)Symmetrizes the plasma in the transverse phase space. For each particle with (x, y, ux, uy), three additional particles are generated with (-x, y, -ux, uy), (x, -y, ux, -uy), and (-x, -y, -ux, -uy). The total number of plasma particles is multiplied by 4. This option is helpful to prevent a numerical seeding of the hosing instability for a plasma with a temperature.
<plasma name> or plasmas.reorder_period
(int) optional (default 0)Reorder particles periodically to speed-up current deposition on GPU for a high-temperature plasma. A good starting point is a period of 4 to reorder plasma particles on every fourth zeta-slice. To disable reordering set this to 0.
<plasma name> or plasmas.n_subcycles
(int) optional (default 1)Number of sub-cycles within the plasma pusher. Currently only implemented for the leapfrog pusher. Must be larger or equal to 1. Sub-cycling is needed if plasma particles move significantly in the transverse direction during a single longitudinal cell. If they move too many cells such that they do not sample certain small transverse structures in the wakefields, sub-cycling is needed and fixes the issue.
<plasma name> or plasmas.reorder_idx_type
(2 int) optional (default 0 0 or 1 1)Change if plasma particles are binned to cells (0), nodes (1) or both (2) for both x and y direction as part of the reordering. The ideal index type depends on the particle shape factor used for deposition. For shape factors 1 and 3, 2^2 and 4^2 cells are deposited per particle respectively, resulting in node centered reordering giving better performance. For shape factors 0 and 2, 1^2 and 3^2 cells are deposited such that cell centered reordering is better. The default is chosen accordingly. If
hipace.depos_derivative_type = 1
, the explicit deposition deposits an additional cell in each direction, making the opposite index type ideal. Since the normal deposition still requires the original index type, the compromise option2 2
can be chosen. This will however require more memory in the binning process.
<plasma name> or plasmas.fine_patch(x,y)
(int) optional (default 0)When using mesh refinement it can be helpful to increase the number of particles per cell drastically in a small part of the domain. For this parameter a function of
x
andy
needs to be specified that evaluates to1
where the number of particles per cell should be higher and0
everywhere else. For example useplasmas.fine_patch(x,y) = "sqrt(x^2+y^2) < 10"
to specify a circle aroundx=0, y=0
with a radius of10
. Note that the function is evaluated at the cell centers of the level zero grid.
<plasma name> or plasmas.fine_ppc
(2 int) optional (default 0 0)The number of plasma particles per cell in x and y inside the fine plasma patch. This must be divisible by the ppc outside the fine patch in both directions.
<plasma name> or plasmas.fine_transition_cells
(int) optional (default 5)Number of cells that are used just outside of the fine plasma patch to smoothly transition between the low and high ppc regions. More transition cells produce less noise but require more particles.
Beam parameters
For the beam parameters, first the names of the beams need to be specified. Afterwards, the beam
parameters for each beam are specified via <beam name>.<beam property> = ...
beams.names
(string) optional (default no_beam)The names of the particle beams, separated by a space. To run without beams, choose the name
no_beam
.
General beam parameters
The general beam parameters are applicable to all particle beam types. More specialized beam parameters,
which are valid only for certain beam types, are introduced further below under
“Option: <injection_type>
”.
<beam name>.injection_type
(string)The injection type for the particle beam. Currently available are
fixed_weight_pdf
,fixed_weight
,fixed_ppc
, andfrom_file
.fixed_weight_pdf
generates a beam with a fixed number of particles with a constant weight where the transverse profile is Gaussian and the longitudinal profile is arbitrary according to a user-specified probability density function. It is more general and faster, and uses less memory thanfixed_weight
.fixed_weight
generates a Gaussian beam with a fixed number of particles with a constant weight.fixed_ppc
generates a beam with a fixed number of particles per cell and varying weights. It can be either a Gaussian or a flattop beam.from_file
reads a beam from openPMD files.
<beam name>.element
(string) optional (default electron)The Physical Element of the plasma. Sets charge, mass and, if available, the specific Ionization Energy of each state. Currently available options are:
electron
,positron
, andproton
.
<beam name>.mass
(float) optional (default m_e)The mass of beam particles. Can also be set with
<beam name>.element
. Must be >0.
<beam name>.charge
(float) optional (default -q_e)The charge of a beam particle. Can also be set with
<beam name>.element
.
<beam name>.n_subcycles
(int) optional (default 10)Number of sub-cycles performed in the beam particle pusher. The particles will be pushed
n_subcycles
times with a time step of dt/n_subcycles. This can be used to improve accuracy in highly non-linear focusing fields.
<beam name> or beams.external_E(x,y,z,t)
(3 float) optional (default 0. 0. 0.)External electric field applied to beam particles as functions of x, y, z and t. The components represent Ex, Ey and Ez respectively. Note that z refers to the location of the beam particle inside the moving frame of reference (zeta) and t to the physical time of the current timestep.
<beam name> or beams.external_B(x,y,z,t)
(3 float) optional (default 0. 0. 0.)External magnetic field applied to beam particles as functions of x, y, z and t. The components represent Bx, By and Bz respectively. Note that z refers to the location of the beam particle inside the moving frame of reference (zeta) and t to the physical time of the current timestep.
<beam name>.do_z_push
(bool) optional (default 1)Whether the beam particles are pushed along the z-axis. The momentum is still fully updated. Note: using
do_z_push = 0
results in unphysical behavior.
<beam name> or beams.do_reset_id_init
(bool) optional (default 0)Whether to reset the ID incrementor to 1 before initializing beam particles.
<beam name> or beams.reorder_period
(int) optional (default 0)Reorder particles periodically to speed-up current deposition and particle push on GPU. A good starting point is a period of 1 to reorder beam particles on every timestep. To disable reordering set this to 0. For very narrow beams the sorting may take longer than the time saved in the beam push and deposition.
<beam name> or beams.reorder_idx_type
(2 int) optional (default 0 0 or 1 1)Change if beam particles are binned to cells (0), nodes (1) or both (2) for both x and y direction as part of the reordering. The ideal index type is different for beam push and beam deposition so some experimentation may be required to find the overall fastest setting for a specific simulation.
Option: fixed_weight_pdf
<beam name>.num_particles
(int)Number of constant weight particles to generate the beam.
<beam name>.pdf
(float)Longitudinal density profile of the beam, given as a probability density function (the transverse profile is Gaussian). This is a parser function of z, giving the charge density integrated in both transverse directions x and y (this is proportional to the beam current profile in the limit \(v_z \simeq c\)). The probability density function is automatically normalized, and combined with
<beam name>.total_charge
or<beam name>.density
within the code to generate the absolute beam profile. Examples (assumingz_center
,z_std
,z_length
,z_slope
,z_min
andz_max
are defined withmy_constants
): - Gaussian:exp(-0.5*((z-z_center)/z_std)^2)
- Cosine:(cos(2*pi*(z-z_center)/z_length)+1)*(2*abs(z-z_center)<z_length)
- Trapezoidal:(z<z_max)*(z>z_min)*(1+z_slope*z)
<beam name>.total_charge
(float)Total charge of the beam (either
total_charge
ordensity
must be specified). Only available when running in SI units. The absolute value of this parameter is used when initializing the beam. Note that in contrast to thefixed_weight
injection type, using<beam name>.radius
or a special pdf to emulatez_min
andz_max
will result in beam particles being redistributed to other locations rather than being deleted. Therefore, the resulting beam will have exactly the specified total charge, but cutting a significant fraction of the charge is not recommended.
<beam name>.density
(float)Peak density of the beam (either
total_charge
ordensity
must be specified). The absolute value of this parameter is used when initializing the beam. Note that this is the peak density of the analytical profile specified by pdf, position_mean and position_std, within the limits of the resolution of the numerical evaluation of the pdf. The actual resulting beam profile consists of randomly distributed particles and will likely feature density fluctuations exceeding the specified peak density.
<beam name>.position_mean
(2 float)The mean position of the beam in
x, y
, separated by a space. Both values can be a function of z. To generate a tilted beam use<beam name>.position_mean = "x_center+(z-z_center)*dx_per_dzeta" "y_center+(z-z_center)*dy_per_dzeta"
.
<beam name>.position_std
(2 float)The rms size of the of the beam in
x, y
, separated by a space. Both values can be a function of z.
<beam name>.u_mean
(3 float)The mean normalized momentum of the beam in
x, y, z
, separated by a space. All values can be a function of z. Normalized momentum is equal to \(= \gamma \beta = \frac{p}{m c}\). An electron beam with a momentum of 1 GeV/c has a u_mean of0 0 1956.951198
while a proton beam with the same momentum has a u_mean of0 0 1.065788933
.
<beam name>.u_std
(3 float)The rms normalized momentum of the beam in
x, y, z
, separated by a space. All values can be a function of z.
<beam name>.do_symmetrize
(bool) optional (default 0)Symmetrizes the beam in the transverse phase space. For each particle with (x, y, ux, uy), three further particles are generated with (-x, y, -ux, uy), (x, -y, ux, -uy), and (-x, -y, -ux, -uy). The total number of particles will still be
beam_name.num_particles
, therefore this option requires that the beam particle number must be divisible by 4.
<beam name>.z_foc
(float) optional (default 0.)Distance at which the beam will be focused, calculated from the position at which the beam is initialized. The beam is assumed to propagate ballistically in-between.
<beam name>.radius
(float) optional (default infinity)Maximum radius
<beam name>.radius
\(= \sqrt{x^2 + y^2}\) within that particles are injected. If<beam name>.density
is specified, beam particles outside of the radius get deleted. If<beam name>.total_charge
is specified, beam particles outside of the radius get new random transverse positions to conserve the total charge.
<beam name>.pdf_ref_ratio
(int) optional (default 4)Into how many segments the pdf is divided per zeta slice for its first-order numerical evaluation.
Option: fixed_weight
<beam name>.num_particles
(int)Number of constant weight particles to generate the beam.
<beam name>.profile
(string) optional (default gaussian)Beam profile. Possible options are
can
(uniform longitudinally, Gaussian transversally) andgaussian
(Gaussian in all directions).
<beam name>.total_charge
(float)Total charge of the beam. Note: Either
total_charge
ordensity
must be specified. The absolute value of this parameter is used when initializing the beam. Note that<beam name>.zmin
,<beam name>.zmax
and<beam name>.radius
can reduce the total charge.
<beam name>.density
(float)Peak density of the beam. Note: Either
total_charge
ordensity
must be specified. The absolute value of this parameter is used when initializing the beam.
<beam name>.position_mean
(3 float)The mean position of the beam in
x, y, z
, separated by a space. The x and y directions can be functions ofz
. To generate a tilted beam use<beam name>.position_mean = "x_center+(z-z_ center)*dx_per_dzeta" "y_center+(z-z_ center)*dy_per_dzeta" "z_center"
.
<beam name>.position_std
(3 float)The rms size of the of the beam in
x, y, z
, separated by a space.
<beam name>.u_mean
(3 float)The mean normalized momentum of the beam in
x, y, z
, separated by a space. Normalized momentum is equal to \(= \gamma \beta = \frac{p}{m c}\). An electron beam with a momentum of 1 GeV/c has a u_mean of0 0 1956.951198
while a proton beam with the same momentum has a u_mean of0 0 1.065788933
.
<beam name>.u_std
(3 float)The rms normalized momentum of the beam in
x, y, z
, separated by a space.
<beam name>.duz_per_uz0_dzeta
(float) optional (default 0.)Relative correlated energy spread per \(\zeta\). Thereby, duz_per_uz0_dzeta * \(\zeta\) * uz_mean is added to uz of the each particle. \(\zeta\) is hereby the particle position relative to the mean longitudinal position of the beam.
<beam name>.do_symmetrize
(bool) optional (default 0)Symmetrizes the beam in the transverse phase space. For each particle with (x, y, ux, uy), three further particles are generated with (-x, y, -ux, uy), (x, -y, ux, -uy), and (-x, -y, -ux, -uy). The total number of particles will still be
beam_name.num_particles
, therefore this option requires that the beam particle number must be divisible by 4.
<beam name>.z_foc
(float) optional (default 0.)Distance at which the beam will be focused, calculated from the position at which the beam is initialized. The beam is assumed to propagate ballistically in-between.
<beam name>.zmin
(float) (default -infinity)Minimum in z at which particles are injected.
<beam name>.zmax
(float) (default infinity)Maximum in z at which particles are injected.
<beam name>.radius
(float) (default infinity)Maximum radius
<beam name>.radius
\(= \sqrt{x^2 + y^2}\) within that particles are injected.
<beam name> or beams.initialize_on_cpu
(bool) optional (default 0)Whether to initialize the beam on the CPU instead of the GPU. Initializing the beam on the CPU can be much slower but is necessary if the full beam does not fit into GPU memory.
Option: fixed_ppc
<beam name>.ppc
(3 int) (default 1 1 1)Number of particles per cell in x-, y-, and z-direction to generate the beam.
<beam name>.profile
(string)Beam profile. Possible options are
flattop
(flat-top radially and longitudinally),gaussian
(Gaussian in all directions), orparsed
(arbitrary analytic function provided by the user). Whenparsed
,<beam name>.density(x,y,z)
must be specified.
<beam name>.density
(float)Peak density of the beam. The absolute value of this parameter is used when initializing the beam.
<beam name>.density(x,y,z)
(float)The density profile of the beam, as a function of spatial dimensions x, y and z. This function uses the parser, see above.
<beam name>.min_density
(float) optional (default 0)Minimum density. Particles with a lower density are not injected. The absolute value of this parameter is used when initializing the beam.
<beam name>.position_mean
(3 float)The mean position of the beam in
x, y, z
, separated by a space.
<beam name>.position_std
(3 float)The rms size of the of the beam in
x, y, z
, separated by a space.
<beam name>.u_mean
(3 float)The mean normalized momentum of the beam in
x, y, z
, separated by a space. Normalized momentum is equal to \(= \gamma \beta = \frac{p}{m c}\). An electron beam with a momentum of 1 GeV/c has a u_mean of0 0 1956.951198
while a proton beam with the same momentum has a u_mean of0 0 1.065788933
.
<beam name>.u_std
(3 float)The rms normalized momentum of the beam in
x, y, z
, separated by a space.
<beam name>.random_ppc
(3 bool) optional (default 0 0 0)Whether the position in (x y z) of the particles is randomized within the cell.
<beam name>.zmin
(float) (default -infinity)Minimum in z at which particles are injected.
<beam name>.zmax
(float) (default infinity)Maximum in z at which particles are injected.
<beam name>.radius
(float) (default infinity)Maximum radius
<beam name>.radius
\(= \sqrt{x^2 + y^2}\) within that particles are injected.
Option: from_file
<beam name> or beams.input_file
(string)Name of the input file. Note: Reading in files with digits in their names (e.g.
openpmd_002135.h5
) can be problematic, it is advised to read them viaopenpmd_%T.h5
and then specify the iteration viabeam_name.iteration = 2135
.
<beam name> or beams.iteration
(integer) optional (default 0)Iteration of the openPMD file to be read in. If the openPMD file contains multiple iterations, or multiple openPMD files are read in, the iteration can be specified. Note: The physical time of the simulation is set to the time of the given iteration (if available).
<beam name>.openPMD_species_name
(string) optional (default <beam name>)Name of the beam to be read in. If an openPMD file contains multiple beams, the name of the beam needs to be specified.
<beam name> or beams.initialize_on_cpu
(bool) optional (default 0)Whether to initialize the beam on the CPU instead of the GPU. Initializing the beam on the CPU can be much slower but is necessary if the full beam does not fit into GPU memory.
SALAME algorithm
HiPACE++ features the Slicing Advanced Loading Algorithm for Minimizing Energy Spread (SALAME) to generate a beam profile that automatically loads the wake optimally, i.e., so that the initial wakefield is flattened by the charge of the beam. Important note: In the algorithm, the weight of the beam particles is adjusted while the plasma response is computed. Since the beam is written to file before the plasma response is calculated, the SALAME beam has incorrect weights in the 0th time step. For more information on the algorithm, see the corresponding publication S. Diederichs et al., Phys. Rev. Accel. Beams 23, 121301 (2020)
<beam name>.do_salame
(bool) optional (default 0)If turned on, the per-slice beam weight in the first time-step is adjusted such that the Ez field is uniform in the beam. This ignores the contributions to jx, jy and rho from the beam in the first time-step. It is recommended to use this option with a fixed weight can beam. If a gaussian beam profile is used, then the zmin and zmax parameters should be used.
hipace.salame_n_iter
(int) optional (default 3)Number of iterations the SALAME algorithm should do when it is used.
hipace.salame_do_advance
(bool) optional (default 1)Whether the SALAME algorithm should calculate the SALAME-beam-only Ez field by advancing plasma (if 1) particles or by approximating it using the chi field (if 0).
hipace.salame_Ez_target(zeta,zeta_initial,Ez_initial)
(string) optional (default Ez_initial)Parser function to specify the target Ez field at the witness beam for SALAME.
zeta
: position of the Ez field to set.zeta_initial
: position where the SALAME algorithm first started.Ez_initial
: field value at zeta_initial. For zeta equal to zeta_initial, the function should return Ez_initial. The default value of this function corresponds to a flat Ez field at the position of the SALAME beam. Note: zeta is always less than or equal to zeta_initial and Ez_initial is typically below zero for electron beams.
Laser parameters
The laser profile is defined by \(a(x,y,z) = a_0 * \mathrm{exp}[-(x^2/w0_x^2 + y^2/w0_y^2 + z^2/L0^2)]\).
The model implemented is the one from [C. Benedetti et al. Plasma Phys. Control. Fusion 60.1: 014002 (2017)].
Unlike for beams
and plasmas
, all the laser pulses are currently stored on the same array,
which you can find in the output openPMD file as a complex array named laserEnvelope.
Parameters starting with lasers.
apply to all laser pulses, parameters starting with <laser name>
apply to a single laser pulse.
lasers.names
(list of string) optional (default no_laser)The names of the laser pulses, separated by a space. To run without a laser, choose the name
no_laser
.
lasers.lambda0
(float)Wavelength of the laser pulses. Currently, all pulses must have the same wavelength.
lasers.use_phase
(bool) optional (default true)Whether the phase terms (\(\theta\) in Eq. (6) of [C. Benedetti et al. Plasma Phys. Control. Fusion 60.1: 014002 (2017)]) are computed and used in the laser envelope advance. Keeping the phase should be more accurate, but can cause numerical issues in the presence of strong depletion/frequency shift.
lasers.solver_type
(string) optional (default multigrid)Type of solver for the laser envelope solver, either
fft
ormultigrid
. Currently, the approximation that the phase is evaluated on-axis only is made with both solvers. With the multigrid solver, we could drop this assumption. For now, the fft solver should be faster, more accurate and more stable, so only use the multigrid one with care.
lasers.MG_tolerance_rel
(float) optional (default 1e-4)Relative error tolerance of the multigrid solver used for the laser pulse.
lasers.MG_tolerance_abs
(float) optional (default 0.)Absolute error tolerance of the multigrid solver used for the laser pulse.
lasers.MG_verbose
(int) optional (default 0)Level of verbosity of the multigrid solver used for the laser pulse.
lasers.MG_average_rhs
(0 or 1) optional (default 1)Whether to use the most stable discretization for the envelope solver.
lasers.input_file
(string) optional (default “”)Path to an openPMD file containing a laser envelope. The file should comply with the LaserEnvelope extension of the openPMD-standard, as generated by LASY. Currently supported geometries: 3D or cylindrical profiles with azimuthal decomposition. The laser pulse is injected in the HiPACE++ simulation so that the beginning of the temporal profile from the file corresponds to the head of the simulation box, and time (in the file) is converted to space (HiPACE++ longitudinal coordinate) with
z = -c*t + const
. If this parameter is set, then the file is used to initialize all lasers instead of using a gaussian profile.
lasers.openPMD_laser_name
(string) optional (default laserEnvelope)Name of the laser envelope field inside the openPMD file to be read in.
lasers.iteration
(int) optional (default 0)Iteration of the openPMD file to be read in.
<laser name>.a0
(float) optional (default 0)Peak normalized vector potential of the laser pulse.
<laser name>.position_mean
(3 float) optional (default 0 0 0)The mean position of the laser in x, y, z.
<laser name>.w0
(2 float) optional (default 0 0)The laser waist in x, y.
<laser name>.L0
(float) optional (default 0)The laser pulse length in z. Use either the pulse length or the pulse duration
<laser name>.tau
.
<laser name>.tau
(float) optional (default 0)The laser pulse duration. The pulse length is set to laser.tau\(*c_0\). Use either the pulse length or the pulse duration.
<laser name>.focal_distance
(float)Distance at which the laser pulse if focused (in the z direction, counted from laser initial position).
<laser name>.propagation_angle_yz
(float)Propagation angle of the pulse in the yz plane (0 is the along the z axis)
Diagnostic parameters
There are different types of diagnostics in HiPACE++. The standard diagnostics are compliant with the openPMD standard. The in-situ diagnostics allow for fast analysis of large beams or the plasma particles.
diagnostic.output_period
(integer) optional (default 0)Output period for standard beam and field diagnostics. Field or beam specific diagnostics can overwrite this parameter. No output is given for
diagnostic.output_period = 0
.
hipace.file_prefix
(string) optional (default diags/hdf5/)Path of the output.
hipace.openpmd_backend
(string) optional (default h5)OpenPMD backend. This can either be
h5
,bp
, orjson
. The default is chosen by what is available. If both Adios2 and HDF5 are available,h5
is used. Note thatjson
is extremely slow and is not recommended for production runs.
Beam diagnostics
diagnostic.beam_output_period
(integer) optional (default 0)Output period for the beam. No output is given for
diagnostic.beam_output_period = 0
. Ifdiagnostic.output_period
is defined, that value is used as the default for this.
diagnostic.beam_data
(string) optional (default all)Names of the beams written to file, separated by a space. The beam names need to be
all
,none
or a subset ofbeams.names
.
Field diagnostics
diagnostic.names
(string) optional (default lev0)The names of all field diagnostics, separated by a space. Multiple diagnostics can be used to limit the output to only a few relevant regions to save on file size. To run without field diagnostics, choose the name
no_field_diag
. If mesh refinement is used, the default becomeslev0 lev1
orlev0 lev1 lev2
.
<diag name> or diagnostic.level
(integer) optional (default 0)From which mesh refinement level the diagnostics should be collected. If
<diag name>
is equal tolev1
, the default for this parameter becomes 1 etc.
<diag name>.output_period
(integer) optional (default 0)Output period for fields. No output is given for
<diag name>.output_period = 0
. Ifdiagnostic.output_period
is defined, that value is used as the default for this.
<diag name> or diagnostic.diag_type
(string)Type of field output. Available options are xyz, xz, yz and xy_integrated. xyz generates a 3D field output. Use 3D output with parsimony, it may increase disk Space usage and simulation time significantly. xz and yz generate 2D field outputs at the center of the y-axis and x-axis, respectively. In case of an even number of grid points, the value is averaged between the two inner grid points. xy_integrated generates 2D field output that has been integrated along the z axis, i.e., it is the sum of the 2D field output over all slices multiplied with dz.
<diag name> or diagnostic.coarsening
(3 int) optional (default 1 1 1)Coarsening ratio of field output in x, y and z direction respectively. The coarsened output is obtained through first order interpolation.
<diag name> or diagnostic.include_ghost_cells
(bool) optional (default 0)Whether the field diagnostics should include ghost cells.
<diag name> or diagnostic.field_data
(string) optional (default all)Names of the fields written to file, separated by a space. The field names need to be
all
,none
or a subset ofExmBy EypBx Ez Bx By Bz Psi
. For the predictor-corrector solver, additionallyjx jy jz rhomjz
are available, which are the current and charge densities of the plasma and the beam, withrhomjz
equal to \(\rho-j_z/c\). For the explicit solver, the current and charge densities of the beam and for all plasmas are separated:jx_beam jy_beam jz_beam
andjx jy rhomjz
are available. Ifrho
is explicitly mentioned asfield_data
, it is deposited by the plasma to be available as a diagnostic. Similarly ifrho_<plasma name>
is explicitly mentioned, the charge density of that plasma species will be separately available as a diagnostic. When a laser pulse is used, the laser complex envelopelaserEnvelope
is available. The plasma proper density (n/gamma) is then also accessible viachi
.
<diag name> or diagnostic.patch_lo
(3 float) optional (default -infinity -infinity -infinity)Lower limit for the diagnostic grid.
<diag name> or diagnostic.patch_hi
(3 float) optional (default infinity infinity infinity)Upper limit for the diagnostic grid.
hipace.deposit_rho
(bool) optional (default 0)If the charge density
rho
of the plasma should be deposited so that it is available as a diagnostic. Otherwise onlyrhomjz
equal to \(\rho-j_z/c\) will be available. Ifrho
is explicitly mentioned indiagnostic.field_data
, then the default will become 1.
hipace.deposit_rho_individual
(bool) optional (default 0)This option works similar to
hipace.deposit_rho
, however the charge density from every plasma species will be deposited into individual fields that are accessible asrho_<plasma name>
indiagnostic.field_data
.
In-situ diagnostics
Besides the standard diagnostics, fast in-situ diagnostics are available. They are most useful when beams with large numbers of particles are used, as the important moments can be calculated in-situ (during the simulation) to largely reduce the simulation’s analysis. In-situ diagnostics compute slice quantities (1 number per quantity per longitudinal cell). For particle beams, they can be used to calculate the main characterizing beam parameters (width, energy spread, emittance, etc.), from which most common beam parameters (e.g. slice and projected emittance, etc.) can be computed. Additionally, the plasma particle properties (e.g, the temperature) can be calculated. For particle quantities, “[…]” stands for averaging over all particles in the current slice; for grid quantities, “[…]” stands for integrating over all cells in the current slice.
For particle beams, the following quantities are calculated per slice and stored:
sum(w), [x], [x^2], [y], [y^2], [z], [z^2], [ux], [ux^2], [uy], [uy^2], [uz], [uz^2], [x*ux], [y*uy], [z*uz], [x*uy], [y*ux], [ux/uz], [uy/uz], [ga], [ga^2], np
.
For plasma particles, the following quantities are calculated per slice and stored:
sum(w), [x], [x^2], [y], [y^2], [ux], [ux^2], [uy], [uy^2], [uz], [uz^2], [ga], [ga^2], np
.
Thereby, “w” stands for weight, “ux” is the normalized momentum in the x direction, “ga” is the Lorentz factor.
Averages and totals over all slices are also provided for convenience under the
respective average
and total
subcategories.
For the field in-situ diagnostics, the following quantities are calculated per slice and stored:
[Ex^2], [Ey^2], [Ez^2], [Bx^2], [By^2], [Bz^2], [ExmBy^2], [EypBx^2], [jz_beam], [Ez*jz_beam]
.
These quantities can be used to calculate the energy stored in the fields.
For the laser in-situ diagnostics, the following quantities are calculated per slice and stored:
max(|a|^2), [|a|^2], [|a|^2*x], [|a|^2*x*x], [|a|^2*y], [|a|^2*y*y], axis(a)
.
Thereby, max(|a|^2)
is the highest value of |a|^2
in the current slice
and axis(a)
gives the complex value of the laser envelope, in the center of every slice.
Additionally, some metadata is also available:
time, step, n_slices, charge, mass, z_lo, z_hi, normalized_density_factor
.
time
and step
refers to the physical time of the simulation and step number of the
current timestep.
n_slices
is the number of slices in the zeta direction.
charge
and mass
relate to a single particle and are for example equal to the
electron charge and mass.
z_lo
and z_hi
are the lower and upper bounds of the z-axis of the simulation domain
specified in the input file and can be used to generate a z/zeta-axis for plotting (note that they corresponds to mesh nodes, while the data is cell-centered).
normalized_density_factor
is equal to dx * dy * dz
in normalized units and 1 in
SI units. It can be used to convert sum(w)
, which specifies the particle density in normalized
units and particle weight in SI units, to the particle weight in both unit systems.
The data is written to a file at <insitu_file_prefix>/reduced_<beam/plasma name>.<MPI rank number>.txt
.
The in-situ diagnostics file format consists of a header part in ASCII containing a JSON object.
When this is parsed into Python it can be converted to a NumPy structured datatype.
The rest of the file, following immediately after the closing }
, is in binary format and
contains all of the in-situ diagnostics along with some metadata. This part can be read using the
structured datatype of the first section.
Use hipace/tools/read_insitu_diagnostics.py
to read the files using this format. Functions to calculate the most useful properties are also provided in that file.
<beam name> or beams.insitu_period
(int) optional (default0
)Period of the beam in-situ diagnostics. 0 means no beam in-situ diagnostics.
<beam name> or beams.insitu_file_prefix
(string) optional (default"diags/insitu"
)Path of the beam in-situ output. Must not be the same as hipace.file_prefix.
<beam name> or beams.insitu_radius
(float) optional (defaultinfinity
)Maximum radius
<beam name>.insitu_radius
\(= \sqrt{x^2 + y^2}\) within which particles are used for the calculation of the insitu diagnostics.
<plasma name> or plasmas.insitu_period
(int) optional (default0
)Period of the plasma in-situ diagnostics. 0 means no plasma in-situ diagnostics.
<plasma name> or plasmas.insitu_file_prefix
(string) optional (default"plasma_diags/insitu"
)Path of the plasma in-situ output. Must not be the same as hipace.file_prefix.
<plasma name> or plasmas.insitu_radius
(float) optional (defaultinfinity
)Maximum radius
<plasma name>.insitu_radius
\(= \sqrt{x^2 + y^2}\) within which particles are used for the calculation of the insitu diagnostics.
fields.insitu_period
(int) optional (default0
)Period of the field in-situ diagnostics. 0 means no field in-situ diagnostics.
fields.insitu_file_prefix
(string) optional (default"diags/field_insitu"
)Path of the field in-situ output. Must not be the same as hipace.file_prefix.
lasers.insitu_period
(int) optional (default0
)Period of the laser in-situ diagnostics. 0 means no laser in-situ diagnostics.
lasers.insitu_file_prefix
(string) optional (default"diags/laser_insitu"
)Path of the laser in-situ output. Must not be the same as hipace.file_prefix.
Additional physics
Additional physics describe the physics modules implemented in HiPACE++ that go beyond the standard electromagnetic equations. This includes ionization (see plasma parameters), binary collisions, and radiation reactions. Since all of these require the actual plasma density, they need a background density in SI units, if the simulation runs in normalized units.
hipace.background_density_SI
(float) optionalBackground plasma density in SI units. Certain physical modules (collisions, ionization, radiation reactions) depend on the actual background density. Hence, in normalized units, they can only be included, if a background plasma density in SI units is provided using this input parameter.
Binary collisions
WARNING: this module is in development.
HiPACE++ proposes an implementation of [Perez et al., Phys. Plasmas 19, 083104 (2012)], inherited from WarpX, for collisions between plasma-plasma and beam-plasma. As collisions depend on the physical density, in normalized units hipace.background_density_SI must be specified.
hipace.collisions
(list of strings) optionalList of names of binary Coulomb collisions. Each will represent collisions between 2 species.
<collision name>.species
(two strings) optionalThe name of the two species for which collisions should be included. This can either be plasma-plasma or beam-plasma collisions. For plasma-plasma collisions, the species can be the same to model collisions within a species. The names must be in plasmas.names or beams.names (for beam-plasma collisions).
<collision name>.CoulombLog
(float) optional (default -1.)Coulomb logarithm used for this collision. If not specified, the Coulomb logarithm is determined from the temperature in each cell.
Radiation reaction
Whether the energy loss due to classical radiation reaction of beam particles is calculated.
<beam name> or beams.do_radiation_reaction
(bool) optional (default 0)Whether the beam particles undergo energy loss due to classical radiation reaction. The implemented radiation reaction model is based on this publication: M. Tamburini et al., NJP 12, 123005 In normalized units, hipace.background_density_SI must be specified.
Get started
If you want to reach out to contributors and users, please consider joining the chat, see Join the community.
See the Build/install HiPACE++ section to install the code locally or on your favorite platform.
A full HiPACE++ input file can be found below, with a detailed description of the main parameters. All the input parameters available in the simulation are available in the Input parameters section.
### AMReX parameters
################################
amr.n_cell = 64 64 100 # number of grid points in x, y, z
amr.max_level = 0 # level of mesh refinement. Currently, only 0 available
#### General simulation settings
################################
max_step = 20 # max time step. 0 computes the fields of the initial beam
diagnostic.output_period = 1 # output period. Last step is always written
hipace.normalized_units = 1 # unit system: SI units: 0, normalized units: 1
hipace.dt = 4.4 # Time step
### Simulation domain
################################
geometry.is_periodic = 1 1 0 # Is periodic?
geometry.prob_lo = -8. -8. -6 # physical domain: dimension must be provided
geometry.prob_hi = 8. 8. 6 # in the respective unit system
### Beam(s)
################################
beams.names = beam # name(s) of the beam(s)
beam.injection_type = fixed_weight # injection type: fixed_weight, fixed_ppc, or from_file
beam.num_particles = 100000 # number of constant weight particles
beam.position_mean = 0. 0. 0 # mean position in x,y,z in the respective unit system
beam.position_std = 0.3 0.3 1.41 # rms size in x,y,z in the respecitve unit system
beam.density = 3. # peak density
beam.u_mean = 0. 0. 2000 # normalized mean momentum in x,y,z
beam.u_std = 0. 0. 0. # normalized rms of the momentum in x,y,z
### Plasma
################################
plasmas.names = plasma # name(s) of the plasma(s)
plasma.element = electron # type of plasma: electron, proton, or an element e.g. He
plasma.density(x,y,z) = 1. # density function in the respective unit systems
plasma.ppc = 1 1 # particles per cell in x,y
### Diagnostics
################################
diagnostic.diag_type = xz # 2D xz slice output. Options: xyz, xz, yz
This input file is inputs_normalized
in examples/get_started/
, which is the recommended place to start running a first simulation.
After compiling HiPACE++ from the HiPACE++ root directory, execute
cd examples/get_started/
../../build/bin/hipace inputs_normalized # run the simulation
Then you can use the openPMD-viewer to read the simulation output data. Directory examples/get_started/ also contains this example notebook
. Reading the simulation output only takes these few Python commands:
# import statements
import numpy as np
import matplotlib.pyplot as plt
from openpmd_viewer import OpenPMDTimeSeries
# Read the simulation data
ts = OpenPMDTimeSeries('./diags/hdf5/')
# Get beam and field data at iteration 20
iteration = 20
x, z = ts.get_particle(species='beam', iteration=iteration, var_list=['x', 'z'])
F, m = ts.get_field(field='Ez', iteration=iteration)
Also please have a look at the production runs for beam-driven and laser-driven wakefield:
Beam-driven wakefield production simulation
max_step = 300
amr.n_cell = 1024 1024 1024
amr.max_level = 0
hipace.max_time = 0.3/clight
diagnostic.output_period = 1
hipace.verbose = 1
hipace.depos_order_xy = 2
hipace.dt = adaptive
hipace.nt_per_betatron = 30
geometry.is_periodic = true true false # Is periodic?
geometry.prob_lo = -250.e-6 -250.e-6 -250.e-6 # physical domain
geometry.prob_hi = 250.e-6 250.e-6 110.e-6
beams.names = driver witness
driver.position_mean = 0. 0. 0.
driver.position_std = 2.e-6 2.e-6 30.e-6
driver.injection_type = fixed_weight
driver.num_particles = 1000000
driver.total_charge = .6e-9
driver.u_mean = 0. 0. 1000.
driver.u_std = 2. 2. 10.
driver.do_symmetrize = 1
witness.position_mean = 0. 0. -160.e-6
witness.position_std = 2.e-6 2.e-6 5.e-6
witness.injection_type = fixed_weight
witness.num_particles = 1000000
witness.total_charge = .2e-9
witness.u_mean = 0. 0. 1000.
witness.u_std = 2. 2. 10.
witness.do_symmetrize = 1
plasmas.names = electron ion
electron.density(x,y,z) = 2.e22
electron.ppc = 1 1
electron.u_mean = 0.0 0.0 0.
electron.element = electron
ion.density(x,y,z) = 2.e22
ion.ppc = 1 1
ion.u_mean = 0.0 0.0 0.
ion.element = H
diagnostic.diag_type = xz
Laser-driven wakefield production simulation
max_step = 500
amr.n_cell = 256 256 1024
my_constants.kp_inv = clight/wp
my_constants.kp = wp/clight
my_constants.wp = sqrt( ne * q_e^2/(m_e *epsilon0) )
my_constants.ne = 1.0505e23
my_constants.Rm = 3*kp_inv
my_constants.Lramp = 6.e-3
hipace.verbose = 1
hipace.dt = 10.*kp_inv/clight
diagnostic.output_period = 10
amr.max_level = 0
geometry.is_periodic = true true false # Is periodic?
geometry.prob_lo = -18*kp_inv -18*kp_inv -7.5*kp_inv
geometry.prob_hi = 18*kp_inv 18*kp_inv 1.5*kp_inv
lasers.names = laser
lasers.lambda0 = 800e-9
lasers.solver_type = multigrid
lasers.MG_tolerance_rel = 1e-5
laser.a0 = 1.9
laser.position_mean = 0. 0. 0
laser.w0 = 3*kp_inv
laser.L0 = 0.5*kp_inv
plasmas.names = plasma
plasma.density(x,y,z) = "ne*(1+4*(x^2+y^2)/(kp^2 * Rm^4 ) ) *
if (z > Lramp, 1, .5*(1-cos(pi*z/Lramp))) *
if (z>0,1,0)"
plasma.ppc = 1 1
plasma.element = electron
plasma.radius = 23.*kp_inv
diagnostic.diag_type = xz
Both use realistic resolutions and run well on modern GPUs (as an estimate, each should take less than 5 min on 16 modern GPUs, e.g., NVIDIA A100).
Additional examples can be found in examples/
.
Data Analysis
Visualize the results
HiPACE++ complies with the openPMD standard, therefore we recommend visualization tools utilizing openPMD.
2D visualization
For 2D visualization, we recommend the openPMD-viewer. Please visit their website for an installation guide and tutorials.
3D visualization
For 3D rendering, we recommend VisualPIC. The latest redesign of VisualPIC is required to visualize openPMD data. Please visit the website for an installation guide.
Community
Join the community
For questions or support, feel free to join the HiPACE++ chat on the Mattermost platform. For your first connection, this will take a few steps:
Click here to reach the login page.
On the “Let’s get started” page, select “Gitlab”. Gitlab is just an auth-provider, you do not need a Gitlab account to join.
On the Helmholtz Codebase page, select Helmholtz AAI at the bottom.
Find your institution and follow the instructions. You should be redirected to a registration form that is pre-filled by your institution. Alternatively, you can select Github and use your Github credentials, provided you have a Github account.
Contribute to HiPACE++
We welcome new contributors!
- To contribute to HiPACE++, the steps are as follows:
Fork the HiPACE++ repo, so you have your own fork
Pull the latest development from baseline, and create your
<new branch>
from itCommit your changes as usual, and push them on your fork
Open a PR between
<new branch>
on your for anddevelopment
on baseline
Documentation
HiPACE++ has a full (functions and classes and their members, albeit sometimes basic) Doxygen-readable documentation. You can compile it with
cd docs
doxygen
open doxyhtml/index.html
The last line would work on MacOS. On another platform, open the html file with your favorite browser.
The HiPACE++ Doxygen documentation can be found here.
Style and conventions
All new element (class, member of a class, struct, function) declared in a .H file must have a Doxygen-readable documentation
Indent four spaces
No tabs allowed
No end-of-line whitespaces allowed
Classes use CamelCase
Objects use snake_case
Lines should not have >100 characters
The declaration and definition of a function should have a space between the function name and the first bracket (
my_function (...)
), function calls should not (my_function(...)
). This is a convention introduce in AMReX sogit grep "my_function ("
returns only the declaration and definition, not the many function calls.
How-to
Make a new release
Find the release tag in all files with something like
git grep '22\.01'
and modify where relevant (be careful with automated search & replace operations, they may cause unwanted changes).- On the main repo page, go to Releases > Draft new release, and
Update the AMReX and openPMD-api versions
Temporarily change the AMReX build version in CMake, similar to the change in file AMReX.cmake in PR 880, so that the HiPACE++ release builds on the corresponding AMReX release.
Click button “Auto-generate release notes” to get a well-formatted list of PRs
Update the commands that you used
Add any additional comments
confirm the release
Revert change in point 2 above so we keep building on AMReX development in-between releases, similar to PR 882.
Once the release is done, Zenodo will generate a DOI. Go to zenodo.org > My Account > Github > HiPACE++, and get the DOI of the last release and copy-paste to the release description