The Simulation class

The Simulation class is top-level class in FBPIC. It contains all the simulation data, and has the high-level method step that performs the PIC cycle.

In addition, its method add_new_species allows to create new particle species, and its method set_moving_window activates the moving window.

class fbpic.main.Simulation(Nz, zmax, Nr, rmax, Nm, dt, p_zmin=-inf, p_zmax=inf, p_rmin=0, p_rmax=inf, p_nz=None, p_nr=None, p_nt=None, n_e=None, zmin=0.0, n_order=-1, dens_func=None, filter_currents=True, v_comoving=None, use_galilean=True, initialize_ions=False, use_cuda=False, n_guard=None, n_damp={'r': 32, 'z': 64}, exchange_period=None, current_correction='curl-free', boundaries={'r': 'reflective', 'z': 'periodic'}, gamma_boost=None, use_all_mpi_ranks=True, particle_shape='linear', verbose_level=1, smoother=None, use_ruyten_shapes=True, use_modified_volume=True)[source]

Initializes a simulation.

By default, this will not create any particle species. You can then add particles species to the simulation by using e.g. the method add_new_species of the simulation object.

Note

As a short-cut, you can also directly create particle species when initializing the Simulation object, by passing the aguments n_e, p_rmin, p_rmax, p_nz, p_nr, p_nt, and dens_func. This will create:

  • an electron species

  • (if initialize_ions is True) an ion species (Hydrogen 1+)

See the docstring of the method add_new_species for the above-mentioned arguments (where n_e has been re-labeled as n).

Parameters:
  • Nz (int) – The number of gridpoints along z

  • Nr (int) – The number of gridpoints along r

  • zmax (float) – The position of the edge of the simulation in z (More precisely, the position of the edge of the last cell)

  • rmax (float) – The position of the edge of the simulation in r (More precisely, the position of the edge of the last cell)

  • Nm (int) – The number of azimuthal modes taken into account. (The simulation uses the modes from m=0 to m=(Nm-1).)

  • dt (float) – The timestep of the simulation

  • n_order (int, optional) – The order of the stencil for z derivatives in the Maxwell solver. Use -1 for infinite order, i.e. for exact dispersion relation in all direction (adviced for single-GPU/single-CPU simulation). Use a positive number (and multiple of 2) for a finite-order stencil (required for multi-GPU/multi-CPU with MPI). A large n_order leads to more overhead in MPI communications, but also to a more accurate dispersion relation for electromagnetic waves. (Typically, n_order = 32 is a good trade-off.) See this article for more information.

  • zmin (float, optional) – The position of the edge of the simulation box. (More precisely, the position of the edge of the first cell)

  • initialize_ions (bool, optional) – Whether to initialize the neutralizing ions

  • filter_currents (bool, optional) – Whether to filter the currents and charge in k space

  • v_comoving (float or None, optional) – If this variable is None, the standard PSATD is used (default). Otherwise, the current is assumed to be “comoving”, i.e. constant with respect to (z - v_comoving * t). This can be done in two ways: either by - Using a PSATD scheme that takes this hypothesis into account - Solving the PSATD scheme in a Galilean frame

  • use_galilean (bool, optional) – Determines which one of the two above schemes is used When use_galilean is true, the whole grid moves with a speed v_comoving

  • use_cuda (bool, optional) – Whether to use CUDA (GPU) acceleration

  • n_guard (int, optional) – Number of guard cells to use at the left and right of a domain, when performing parallel (MPI) computation or when using open boundaries. Defaults to None, which calculates the required guard cells for n_order automatically (approx 2*n_order). If no MPI is used and in the case of open boundaries with an infinite order stencil, n_guard defaults to 64, if not set otherwise.

  • n_damp (dict, optional) – A dictionary with ‘z’ and ‘r’ as keys, and integers as values. The integers represent the number of damping cells in the longitudinal (z) and transverse (r) directions, respectively. The damping cells in z are only used if boundaries[‘z’] is ‘open’, and are added at the left and right edge of the simulation domain. The damping cells in r are used only if boundaries[‘r’] is ‘open’, and are added at upper radial boundary (at rmax).

  • exchange_period (int, optional) – Number of iterations before which the particles are exchanged. If set to None, the maximum exchange period is calculated automatically: Within exchange_period timesteps, the particles should never be able to travel more than (n_guard/2 - particle_shape order) cells. (Setting exchange_period to small values can substantially affect the performance)

  • boundaries (dict, optional) –

    A dictionary with ‘z’ and ‘r’ as keys, and strings as values. This specifies the field boundary in the longitudinal (z) and transverse (r) direction respectively:

    • boundaries[‘z’] can be either ‘periodic’ or ‘open’ (for field-absorbing boundary).

    • boundaries[‘r’] can be either ‘reflective’ or ‘open’ (for field-absorbing boundary). For ‘open’, this adds Perfectly-Matched-Layers in the radial direction ; note that the computation is significantly more costly in this case.

  • current_correction (string, optional) – The method used in order to ensure that the continuity equation is satisfied. Either curl-free or cross-deposition. curl-free is faster but less local.

  • gamma_boost (float, optional) – When running the simulation in a boosted frame, set the value of gamma_boost to the corresponding Lorentz factor. All the other quantities (zmin, zmax, n_e, etc.) are to be given in the lab frame.

  • use_all_mpi_ranks (bool, optional) –

    When launching the simulation with mpirun:

    • if use_all_mpi_ranks is True (default): All the MPI ranks will contribute to the same simulation, using domain-decomposition to share the work.

    • if use_all_mpi_ranks is False: Each MPI rank will run an independent simulation. This can be useful when running parameter scans. In this case, make sure that your input script is written so that the input parameters and output folder depend on the MPI rank.

  • particle_shape (str, optional) – Set the particle shape for the charge/current deposition. Possible values are ‘cubic’, ‘linear’. (‘cubic’ corresponds to third order shapes and ‘linear’ to first order shapes).

  • verbose_level (int, optional) – Print information about the simulation setup after initialization of the Simulation class. 0 - Print no information 1 (Default) - Print basic information 2 - Print detailed information

  • smoother (an instance of BinomialSmoother, optional) – Determines how the charge and currents are smoothed. (Default: one-pass binomial filter and no compensator.)

  • use_ruyten_shapes (bool, optional) – Whether to use Ruyten shape factors for the particle deposition. (Ruyten JCP 105 (1993) https://doi.org/10.1006/jcph.1993.1070) This ensures that a uniform distribution of macroparticles leads to a uniform charge density deposited on the grid, even close to the axis (in the limit of many particles in r).

  • use_modified_volume (bool, optional) – Whether to use a slightly-modified, effective cell volume, that ensures that the charge deposited near the axis is correctly taken into account by the spectral cylindrical Maxwell solver.

add_new_species(q, m, n=None, dens_func=None, p_nz=None, p_nr=None, p_nt=None, p_zmin=-inf, p_zmax=inf, p_rmin=0, p_rmax=inf, uz_m=0.0, ux_m=0.0, uy_m=0.0, uz_th=0.0, ux_th=0.0, uy_th=0.0, continuous_injection=True, boost_positions_in_dens_func=False, is_tracer=False)[source]

Create a new species (i.e. an instance of Particles) with charge q and mass m. Add it to the simulation (i.e. to the list Simulation.ptcl), and return it, so that the methods of the Particles class can be used, for this particular species.

In addition, if n is set, then new macroparticles will be created within this species (in an evenly-spaced manner).

For boosted-frame simulations (i.e. where gamma_boost as been passed to the Simulation object), all quantities that are explicitly mentioned to be in the lab frame below are automatically converted to the boosted frame.

Note

For the arguments below, it is recommended to have at least p_nt = 4*Nm (except in the case Nm=1, for which p_nt=1 is sufficient). In other words, the required number of macroparticles along theta (in order for the simulation to be properly resolved) increases with the number of azimuthal modes used.

Parameters:
  • q (float (in Coulombs)) – Charge of the particle species

  • m (float (in kg)) – Mass of the particle species

  • n (float (in particles per m^3) or None, optional) – Density of physical particles (in the lab frame). If this is None, no macroparticles will be created. If n is not None, evenly-spaced macroparticles will be generated.

  • dens_func (callable, optional) –

    A function of the form dens_func( z, r ) where z and r are 1d arrays, or dens( x, y, z) where x, y and z are 1d arrays, and which returns a 1d array containing the density relative to `n` (i.e. a number between 0 and 1) at the given positions

    For boosted-frame simulation: if you set boost_positions_in_dens_func to True, then dens_func can be expressed as a function of z taken in the lab frame. Otherwise, it has to be expressed as a function of z taken in the boosted frame.

  • p_nz (int, optional) – The number of macroparticles per cell along the z direction

  • p_nr (int, optional) – The number of macroparticles per cell along the r direction

  • p_nt (int, optional) – The number of macroparticles along the theta direction

  • p_zmin (float (in meters), optional) – The minimal z position above which the particles are initialized (in the lab frame). Default: left edge of the simulation box.

  • p_zmax (float (in meters), optional) – The maximal z position below which the particles are initialized (in the lab frame). Default: right edge of the simulation box.

  • p_rmin (float (in meters), optional) – The minimal r position above which the particles are initialized (in the lab frame). Default: 0

  • p_rmax (floats (in meters), optional) – The maximal r position below which the particles are initialized (in the lab frame). Default: upper edge of the simulation box.

  • uz_m (floats (dimensionless), optional) – Normalized mean momenta (in the lab frame) of the injected particles in each direction

  • ux_m (floats (dimensionless), optional) – Normalized mean momenta (in the lab frame) of the injected particles in each direction

  • uy_m (floats (dimensionless), optional) – Normalized mean momenta (in the lab frame) of the injected particles in each direction

  • uz_th (floats (dimensionless), optional) – Normalized thermal momenta (in the lab frame) in each direction

  • ux_th (floats (dimensionless), optional) – Normalized thermal momenta (in the lab frame) in each direction

  • uy_th (floats (dimensionless), optional) – Normalized thermal momenta (in the lab frame) in each direction

  • continuous_injection (bool, optional) – Whether to continuously inject the particles, in the case of a moving window

  • boost_positions_in_dens_func (bool, optional) – For boosted-frame simulations: whether to automatically take into account the Lorentz transformation of the positions, in dens_func

  • is_tracer (bool, optional) – Setting this flag to True will allow this particle to move as a normal particle with given mass and charge, but will generate no current. This allows them to be passive tracers inside the plasma.

Returns:

new_species

Return type:

an instance of the Particles class

set_moving_window(v=299792458.0, ux_m=None, uy_m=None, uz_m=None, ux_th=None, uy_th=None, uz_th=None, gamma_boost=None)[source]

Initializes a moving window for the simulation.

Parameters:
  • v (float (in meters per seconds), optional) – The speed of the moving window

  • ux_m (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • uy_m (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • uz_m (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • ux_th (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • uy_th (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • uz_th (float (dimensionless), optional) – Unused, kept for backward-compatibility

  • gamma_boost (float, optional) – Unused; kept for backward compatibility

step(N=1, correct_currents=True, correct_divE=False, use_true_rho=False, move_positions=True, move_momenta=True, show_progress=True)[source]

Perform N PIC cycles.

Parameters:
  • N (int, optional) – The number of timesteps to take Default: N=1

  • correct_currents (bool, optional) – Whether to correct the currents in spectral space

  • correct_divE (bool, optional) – Whether to correct the divergence of E in spectral space

  • use_true_rho (bool, optional) – Whether to use the true rho deposited on the grid for the field push or not. (requires initialize_ions = True)

  • move_positions (bool, optional) – Whether to move or freeze the particles’ positions

  • move_momenta (bool, optional) – Whether to move or freeze the particles’ momenta

  • show_progress (bool, optional) – Whether to show a progression bar

fbpic.utils.random_seed.set_random_seed(random_seed)[source]

Fix the seed of the random number generators.

Fixing a seed helps ensure that repeatedly running the same simulation gives the same result (despite the Monte Carlo parts of the code, e.g. ionization, gaussian beam generation, etc.)

random_seed: int

The seed of the random number generator.