Python (PICMI)

WarpX uses the PICMI standard for its Python input files. Python version 3.7 or newer is required.

Example input files can be found in the examples section. In the input file, instances of classes are created defining the various aspects of the simulation. The Simulation object is the central object, where the instances are passed, defining the simulation time, field solver, registered species, etc.

Classes

Simulation and grid setup

Simulation

class pywarpx.picmi.Simulation(solver=None, time_step_size=None, max_steps=None, max_time=None, verbose=None, particle_shape='linear', gamma_boost=None, cpu_split=None, load_balancing=None, **kw)[source]

Creates a Simulation object

Parameters
  • solver (field solver instance) – This is the field solver to be used in the simulation. It should be an instance of field solver classes.

  • time_step_size (float) – Absolute time step size of the simulation [s]. Needed if the CFL is not specified elsewhere.

  • max_steps (integer) – Maximum number of time steps. Specify either this, or max_time, or use the step function directly.

  • max_time (float) – Maximum physical time to run the simulation [s]. Specify either this, or max_steps, or use the step function directly.

  • verbose (integer, optional) – Verbosity flag. A larger integer results in more verbose output

  • particle_shape ({'NGP', 'linear', 'quadratic', 'cubic'}) – Default particle shape for species added to this simulation

  • gamma_boost (float, optional) – Lorentz factor of the boosted simulation frame. Note that all input values should be in the lab frame.

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_current_deposition_algo ({'direct', 'esirkepov', and 'vay'}, optional) – Current deposition algorithm. The default depends on conditions.

  • warpx_charge_deposition_algo ({'standard'}, optional) – Charge deposition algorithm.

  • warpx_field_gathering_algo ({'energy-conserving', 'momentum-conserving'}, optional) – Field gathering algorithm. The default depends on conditions.

  • warpx_particle_pusher_algo ({'boris', 'vay', 'higuera'}, default='boris') – Particle pushing algorithm.

  • warpx_use_filter (bool, optional) – Whether to use filtering. The default depends on the conditions.

  • warpx_serialize_initial_conditions (bool, default=False) – Controls the random numbers used for initialization. This parameter should only be used for testing and continuous integration.

  • warpx_do_dynamic_scheduling (bool, default=True) – Whether to do dynamic scheduling with OpenMP

  • warpx_load_balance_intervals (string, default='0') – The intervals for doing load balancing

  • warpx_load_balance_efficiency_ratio_threshold (float, default=1.1) – (See documentation)

  • warpx_load_balance_with_sfc (bool, default=0) – (See documentation)

  • warpx_load_balance_knapsack_factor (float, default=1.24) – (See documentation)

  • warpx_load_balance_costs_update ({'heuristic' or 'timers' or 'gpuclock'}, optional) – (See documentation)

  • warpx_costs_heuristic_particles_wt (float, optional) – (See documentation)

  • warpx_costs_heuristic_cells_wt (float, optional) – (See documentation)

  • warpx_use_fdtd_nci_corr (bool, optional) – Whether to use the NCI correction when using the FDTD solver

  • warpx_amr_check_input (bool, optional) – Whether AMReX should perform checks on the input (primarily related to the max grid size and blocking factors)

  • warpx_amr_restart (string, optional) – The name of the restart to use

  • warpx_zmax_plasma_to_compute_max_step (float, optional) – Sets the simulation run time based on the maximum z value

  • warpx_collisions (collision instance, optional) – The collision instance specifying the particle collisions

  • warpx_embedded_boundary (embedded boundary instance, optional) –

  • warpx_break_signals (list of strings) – Signals on which to break

  • warpx_checkpoint_signals (list of strings) – Signals on which to write out a checkpoint

add_laser(laser, injection_method)

Add a laser pulses that to be injected in the simulation

Parameters
  • laser_profile (laser instance) – One of laser profile instances. Specifies the physical properties of the laser pulse (e.g. spatial and temporal profile, wavelength, amplitude, etc.).

  • injection_method (laser injection instance, optional) – Specifies how the laser is injected (numerically) into the simulation (e.g. through a laser antenna, or directly added to the mesh). This argument describes an algorithm, not a physical object. It is up to each code to define the default method of injection, if the user does not provide injection_method.

add_species(species, layout, initialize_self_field=None)

Add species to be used in the simulation

Parameters
  • species (species instance) – An instance of one of the PICMI species objects. Defines species to be added from the physical point of view (e.g. charge, mass, initial distribution of particles).

  • layout (layout instance) – An instance of one of the PICMI particle layout objects. Defines how particles are added into the simulation, from the numerical point of view.

  • initialize_self_field (bool, optional) – Whether the initial space-charge fields of this species is calculated and added to the simulation

step(nsteps=None, mpi_comm=None)[source]

Run the simulation for nsteps timesteps

Parameters

nsteps (integer, default=1) – The number of timesteps

write_input_file(file_name='inputs')[source]

Write the parameters of the simulation, as defined in the PICMI input, into a code-specific input file.

This can be used for codes that are not Python-driven (e.g. compiled, pure C++ or Fortran codes) and expect a text input in a given format.

Parameters

file_name (string) – The path to the file that will be created

Constants

For convenience, the PICMI interface defines the following constants, which can be used directly inside any PICMI script. The values are in SI units.

  • picmi.constants.c: The speed of light in vacuum.

  • picmi.constants.ep0: The vacuum permittivity \(\epsilon_0\)

  • picmi.constants.mu0: The vacuum permeability \(\mu_0\)

  • picmi.constants.q_e: The elementary charge (absolute value of the charge of an electron).

  • picmi.constants.m_e: The electron mass

  • picmi.constants.m_p: The proton mass

Field solvers define the updates of electric and magnetic fields.

ElectromagneticSolver

class pywarpx.picmi.ElectromagneticSolver(grid, method=None, stencil_order=None, cfl=None, l_nodal=None, source_smoother=None, field_smoother=None, subcycling=None, galilean_velocity=None, divE_cleaning=None, divB_cleaning=None, pml_divE_cleaning=None, pml_divB_cleaning=None, **kw)[source]

Electromagnetic field solver

Parameters
  • grid (grid instance) – Grid object for the diagnostic

  • method ({'Yee', 'CKC', 'Lehe', 'PSTD', 'PSATD', 'GPSTD', 'DS', 'ECT'}) –

    The advance method use to solve Maxwell’s equations. The default method is code dependent.

  • stencil_order (vector of integers) – Order of stencil for each axis (-1=infinite)

  • cfl (float, optional) – Fraction of the Courant-Friedrich-Lewy criteria [1]

  • l_nodal (bool, optional) – Quantities are at nodes if True, staggered otherwise

  • source_smoother (smoother instance, optional) – Smoother object to apply to the sources

  • field_smoother (smoother instance, optional) – Smoother object to apply to the fields

  • subcycling (integer, optional) – Level of subcycling for the GPSTD solver

  • galilean_velocity (vector of floats, optional) – Velocity of Galilean reference frame [m/s]

  • divE_cleaning (bool, optional) – Solver uses div(E) cleaning if True

  • divB_cleaning (bool, optional) – Solver uses div(B) cleaning if True

  • pml_divE_cleaning (bool, optional) – Solver uses div(E) cleaning in the PML if True

  • pml_divB_cleaning (bool, optional) – Solver uses div(B) cleaning in the PML if True

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_pml_ncell (integer, optional) – The depth of the PML, in number of cells

  • warpx_periodic_single_box_fft (bool, default=False) – Whether to do the spectral solver FFTs assuming a single simulation block

  • warpx_current_correction (bool, default=True) – Whether to do the current correction for the spectral solver. See documentation for exceptions to the default value.

  • warpx_psatd_update_with_rho (bool, optional) – Whether to update with the actual rho for the spectral solver

  • warpx_psatd_do_time_averaging (bool, optional) – Whether to do the time averaging for the spectral solver

  • warpx_do_pml_in_domain (bool, default=False) – Whether to do the PML boundaries within the domain (versus in the guard cells)

  • warpx_pml_has_particles (bool, default=False) – Whether to allow particles in the PML region

  • warpx_do_pml_j_damping (bool, default=False) – Whether to do damping of J in the PML

ElectrostaticSolver

class pywarpx.picmi.ElectrostaticSolver(grid, method=None, required_precision=None, maximum_iterations=None, **kw)[source]

Electrostatic field solver

Parameters
  • grid (grid instance) – Grid object for the diagnostic

  • method (string) – One of ‘FFT’, or ‘Multigrid’

  • required_precision (float, optional) – Level of precision required for iterative solvers

  • maximum_iterations (integer, optional) – Maximum number of iterations for iterative solvers

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_relativistic (bool, default=False) – Whether to use the relativistic solver or lab frame solver

  • warpx_absolute_tolerance (float, default=0.) – Absolute tolerance on the lab fram solver

  • warpx_self_fields_verbosity (integer, default=2) – Level of verbosity for the lab frame solver

Cartesian3DGrid

class pywarpx.picmi.Cartesian3DGrid(number_of_cells=None, lower_bound=None, upper_bound=None, lower_boundary_conditions=None, upper_boundary_conditions=None, nx=None, ny=None, nz=None, xmin=None, xmax=None, ymin=None, ymax=None, zmin=None, zmax=None, bc_xmin=None, bc_xmax=None, bc_ymin=None, bc_ymax=None, bc_zmin=None, bc_zmax=None, moving_window_velocity=None, refined_regions=[], lower_bound_particles=None, upper_bound_particles=None, xmin_particles=None, xmax_particles=None, ymin_particles=None, ymax_particles=None, zmin_particles=None, zmax_particles=None, lower_boundary_conditions_particles=None, upper_boundary_conditions_particles=None, bc_xmin_particles=None, bc_xmax_particles=None, bc_ymin_particles=None, bc_ymax_particles=None, bc_zmin_particles=None, bc_zmax_particles=None, guard_cells=None, pml_cells=None, **kw)[source]

Three dimensional Cartesian grid Parameters can be specified either as vectors or separately. (If both are specified, the vector is used.)

Parameters
  • number_of_cells (vector of integers) – Number of cells along each axis (number of nodes is number_of_cells+1)

  • lower_bound (vector of floats) – Position of the node at the lower bound [m]

  • upper_bound (vector of floats) – Position of the node at the upper bound [m]

  • lower_boundary_conditions (vector of strings) – Conditions at lower boundaries, periodic, open, dirichlet, or neumann

  • upper_boundary_conditions (vector of strings) – Conditions at upper boundaries, periodic, open, dirichlet, or neumann

  • nx (integer) – Number of cells along X (number of nodes=nx+1)

  • ny (integer) – Number of cells along Y (number of nodes=ny+1)

  • nz (integer) – Number of cells along Z (number of nodes=nz+1)

  • xmin (float) – Position of first node along X [m]

  • xmax (float) – Position of last node along X [m]

  • ymin (float) – Position of first node along Y [m]

  • ymax (float) – Position of last node along Y [m]

  • zmin (float) – Position of first node along Z [m]

  • zmax (float) – Position of last node along Z [m]

  • bc_xmin (string) – Boundary condition at min X: One of periodic, open, dirichlet, or neumann

  • bc_xmax (string) – Boundary condition at max X: One of periodic, open, dirichlet, or neumann

  • bc_ymin (string) – Boundary condition at min Y: One of periodic, open, dirichlet, or neumann

  • bc_ymax (string) – Boundary condition at max Y: One of periodic, open, dirichlet, or neumann

  • bc_zmin (string) – Boundary condition at min Z: One of periodic, open, dirichlet, or neumann

  • bc_zmax (string) – Boundary condition at max Z: One of periodic, open, dirichlet, or neumann

  • moving_window_velocity (vector of floats, optional) – Moving frame velocity [m/s]

  • refined_regions (list of lists, optional) – List of refined regions, each element being a list of the format [level, lo, hi, refinement_factor], with level being the refinement level, with 1 being the first level of refinement, 2 being the second etc, lo and hi being vectors of length 3 specifying the extent of the region, and refinement_factor defaulting to [2,2,2] (relative to next lower level)

  • lower_bound_particles (vector of floats, optional) – Position of particle lower bound [m]

  • upper_bound_particles (vector of floats, optional) – Position of particle upper bound [m]

  • xmin_particles (float, optional) – Position of min particle boundary along X [m]

  • xmax_particles (float, optional) – Position of max particle boundary along X [m]

  • ymin_particles (float, optional) – Position of min particle boundary along Y [m]

  • ymax_particles (float, optional) – Position of max particle boundary along Y [m]

  • float (zmin_particles) – Position of min particle boundary along Z [m]

  • optional – Position of min particle boundary along Z [m]

  • zmax_particles (float, optional) – Position of max particle boundary along Z [m]

  • lower_boundary_conditions_particles (vector of strings, optional) – Conditions at lower boundaries for particles, periodic, absorbing, reflect or thermal

  • upper_boundary_conditions_particles (vector of strings, optional) – Conditions at upper boundaries for particles, periodic, absorbing, reflect or thermal

  • bc_xmin_particles (string, optional) – Boundary condition at min X for particles: One of periodic, absorbing, reflect, thermal

  • bc_xmax_particles (string, optional) – Boundary condition at max X for particles: One of periodic, absorbing, reflect, thermal

  • bc_ymin_particles (string, optional) – Boundary condition at min Y for particles: One of periodic, absorbing, reflect, thermal

  • bc_ymax_particles (string, optional) – Boundary condition at max Y for particles: One of periodic, absorbing, reflect, thermal

  • bc_zmin_particles (string, optional) – Boundary condition at min Z for particles: One of periodic, absorbing, reflect, thermal

  • bc_zmax_particles (string, optional) – Boundary condition at max Z for particles: One of periodic, absorbing, reflect, thermal

  • guard_cells (vector of integers, optional) – Number of guard cells used along each direction

  • pml_cells (vector of integers, optional) – Number of Perfectly Matched Layer (PML) cells along each direction

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_max_grid_size (integer, default=32) – Maximum block size in either direction

  • warpx_max_grid_size_x (integer, optional) – Maximum block size in x direction

  • warpx_max_grid_size_y (integer, optional) – Maximum block size in z direction

  • warpx_max_grid_size_z (integer, optional) – Maximum block size in z direction

  • warpx_blocking_factor (integer, optional) – Blocking factor (which controls the block size)

  • warpx_blocking_factor_x (integer, optional) – Blocking factor (which controls the block size) in the x direction

  • warpx_blocking_factor_y (integer, optional) – Blocking factor (which controls the block size) in the z direction

  • warpx_blocking_factor_z (integer, optional) – Blocking factor (which controls the block size) in the z direction

  • warpx_potential_lo_x (float, default=0.) – Electrostatic potential on the lower x boundary

  • warpx_potential_hi_x (float, default=0.) – Electrostatic potential on the upper x boundary

  • warpx_potential_lo_y (float, default=0.) – Electrostatic potential on the lower z boundary

  • warpx_potential_hi_y (float, default=0.) – Electrostatic potential on the upper z boundary

  • warpx_potential_lo_z (float, default=0.) – Electrostatic potential on the lower z boundary

  • warpx_potential_hi_z (float, default=0.) – Electrostatic potential on the upper z boundary

Cartesian2DGrid

class pywarpx.picmi.Cartesian2DGrid(number_of_cells=None, lower_bound=None, upper_bound=None, lower_boundary_conditions=None, upper_boundary_conditions=None, nx=None, ny=None, xmin=None, xmax=None, ymin=None, ymax=None, bc_xmin=None, bc_xmax=None, bc_ymin=None, bc_ymax=None, moving_window_velocity=None, refined_regions=[], lower_bound_particles=None, upper_bound_particles=None, xmin_particles=None, xmax_particles=None, ymin_particles=None, ymax_particles=None, lower_boundary_conditions_particles=None, upper_boundary_conditions_particles=None, bc_xmin_particles=None, bc_xmax_particles=None, bc_ymin_particles=None, bc_ymax_particles=None, guard_cells=None, pml_cells=None, **kw)[source]

Two dimensional Cartesian grid Parameters can be specified either as vectors or separately. (If both are specified, the vector is used.)

Parameters
  • number_of_cells (vector of integers) – Number of cells along each axis (number of nodes is number_of_cells+1)

  • lower_bound (vector of floats) – Position of the node at the lower bound [m]

  • upper_bound (vector of floats) – Position of the node at the upper bound [m]

  • lower_boundary_conditions (vector of strings) – Conditions at lower boundaries, periodic, open, dirichlet, or neumann

  • upper_boundary_conditions (vector of strings) – Conditions at upper boundaries, periodic, open, dirichlet, or neumann

  • nx (integer) – Number of cells along X (number of nodes=nx+1)

  • ny (integer) – Number of cells along Y (number of nodes=ny+1)

  • xmin (float) – Position of first node along X [m]

  • xmax (float) – Position of last node along X [m]

  • ymin (float) – Position of first node along Y [m]

  • ymax (float) – Position of last node along Y [m]

  • bc_xmin (vector of strings) – Boundary condition at min X: One of periodic, open, dirichlet, or neumann

  • bc_xmax (vector of strings) – Boundary condition at max X: One of periodic, open, dirichlet, or neumann

  • bc_ymin (vector of strings) – Boundary condition at min Y: One of periodic, open, dirichlet, or neumann

  • bc_ymax (vector of strings) – Boundary condition at max Y: One of periodic, open, dirichlet, or neumann

  • moving_window_velocity (vector of floats, optional) – Moving frame velocity [m/s]

  • refined_regions (list of lists, optional) – List of refined regions, each element being a list of the format [level, lo, hi, refinement_factor], with level being the refinement level, with 1 being the first level of refinement, 2 being the second etc, lo and hi being vectors of length 2 specifying the extent of the region, and refinement_factor defaulting to [2,2] (relative to next lower level)

  • lower_bound_particles (vector of floats, optional) – Position of particle lower bound [m]

  • upper_bound_particles (vector of floats, optional) – Position of particle upper bound [m]

  • xmin_particles (float, optional) – Position of min particle boundary along X [m]

  • xmax_particles (float, optional) – Position of max particle boundary along X [m]

  • ymin_particles (float, optional) – Position of min particle boundary along Y [m]

  • ymax_particles (float, optional) – Position of max particle boundary along Y [m]

  • lower_boundary_conditions_particles (vector of strings, optional) – Conditions at lower boundaries for particles, periodic, absorbing, reflect or thermal

  • upper_boundary_conditions_particles (vector of strings, optional) – Conditions at upper boundaries for particles, periodic, absorbing, reflect or thermal

  • bc_xmin_particles (string, optional) – Boundary condition at min X for particles: One of periodic, absorbing, reflect, thermal

  • bc_xmax_particles (string, optional) – Boundary condition at max X for particles: One of periodic, absorbing, reflect, thermal

  • bc_ymin_particles (string, optional) – Boundary condition at min Y for particles: One of periodic, absorbing, reflect, thermal

  • bc_ymax_particles (string, optional) – Boundary condition at max Y for particles: One of periodic, absorbing, reflect, thermal

  • guard_cells (vector of integers, optional) – Number of guard cells used along each direction

  • pml_cells (vector of integers, optional) – Number of Perfectly Matched Layer (PML) cells along each direction

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_max_grid_size (integer, default=32) – Maximum block size in either direction

  • warpx_max_grid_size_x (integer, optional) – Maximum block size in x direction

  • warpx_max_grid_size_y (integer, optional) – Maximum block size in z direction

  • warpx_blocking_factor (integer, optional) – Blocking factor (which controls the block size)

  • warpx_blocking_factor_x (integer, optional) – Blocking factor (which controls the block size) in the x direction

  • warpx_blocking_factor_y (integer, optional) – Blocking factor (which controls the block size) in the z direction

  • warpx_potential_lo_x (float, default=0.) – Electrostatic potential on the lower x boundary

  • warpx_potential_hi_x (float, default=0.) – Electrostatic potential on the upper x boundary

  • warpx_potential_lo_z (float, default=0.) – Electrostatic potential on the lower z boundary

  • warpx_potential_hi_z (float, default=0.) – Electrostatic potential on the upper z boundary

Cartesian1DGrid

class pywarpx.picmi.Cartesian1DGrid(number_of_cells=None, lower_bound=None, upper_bound=None, lower_boundary_conditions=None, upper_boundary_conditions=None, nx=None, xmin=None, xmax=None, bc_xmin=None, bc_xmax=None, moving_window_velocity=None, refined_regions=[], lower_bound_particles=None, upper_bound_particles=None, xmin_particles=None, xmax_particles=None, lower_boundary_conditions_particles=None, upper_boundary_conditions_particles=None, bc_xmin_particles=None, bc_xmax_particles=None, guard_cells=None, pml_cells=None, **kw)[source]

One-dimensional Cartesian grid Parameters can be specified either as vectors or separately. (If both are specified, the vector is used.)

Parameters
  • number_of_cells (vector of integers) – Number of cells along each axis (number of nodes is number_of_cells+1)

  • lower_bound (vector of floats) – Position of the node at the lower bound [m]

  • upper_bound (vector of floats) – Position of the node at the upper bound [m]

  • lower_boundary_conditions (vector of strings) – Conditions at lower boundaries, periodic, open, dirichlet, or neumann

  • upper_boundary_conditions (vector of strings) – Conditions at upper boundaries, periodic, open, dirichlet, or neumann

  • nx (integer) – Number of cells along X (number of nodes=nx+1)

  • xmin (float) – Position of first node along X [m]

  • xmax (float) – Position of last node along X [m]

  • bc_xmin (vector of strings) – Boundary condition at min X: One of periodic, open, dirichlet, or neumann

  • bc_xmax (vector of strings) – Boundary condition at max X: One of periodic, open, dirichlet, or neumann

  • moving_window_velocity (vector of floats, optional) – Moving frame velocity [m/s]

  • refined_regions (list of lists, optional) – List of refined regions, each element being a list of the format [level, lo, hi, refinement_factor], with level being the refinement level, with 1 being the first level of refinement, 2 being the second etc, lo and hi being vectors of length 2 specifying the extent of the region, and refinement_factor defaulting to [2,2] (relative to next lower level)

  • lower_bound_particles (vector of floats, optional) – Position of particle lower bound [m]

  • upper_bound_particles (vector of floats, optional) – Position of particle upper bound [m]

  • xmin_particles (float, optional) – Position of min particle boundary along X [m]

  • xmax_particles (float, optional) – Position of max particle boundary along X [m]

  • lower_boundary_conditions_particles (vector of strings, optional) – Conditions at lower boundaries for particles, periodic, absorbing, reflect or thermal

  • upper_boundary_conditions_particles (vector of strings, optional) – Conditions at upper boundaries for particles, periodic, absorbing, reflect or thermal

  • bc_xmin_particles (string, optional) – Boundary condition at min X for particles: One of periodic, absorbing, reflect, thermal

  • bc_xmax_particles (string, optional) – Boundary condition at max X for particles: One of periodic, absorbing, reflect, thermal

  • guard_cells (vector of integers, optional) – Number of guard cells used along each direction

  • pml_cells (vector of integers, optional) – Number of Perfectly Matched Layer (PML) cells along each direction

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_max_grid_size (integer, default=32) – Maximum block size in either direction

  • warpx_max_grid_size_x (integer, optional) – Maximum block size in longitudinal direction

  • warpx_blocking_factor (integer, optional) – Blocking factor (which controls the block size)

  • warpx_blocking_factor_x (integer, optional) – Blocking factor (which controls the block size) in the longitudinal direction

  • warpx_potential_lo_z (float, default=0.) – Electrostatic potential on the lower longitudinal boundary

  • warpx_potential_hi_z (float, default=0.) – Electrostatic potential on the upper longitudinal boundary

CylindricalGrid

class pywarpx.picmi.CylindricalGrid(number_of_cells=None, lower_bound=None, upper_bound=None, lower_boundary_conditions=None, upper_boundary_conditions=None, nr=None, nz=None, n_azimuthal_modes=None, rmin=None, rmax=None, zmin=None, zmax=None, bc_rmin=None, bc_rmax=None, bc_zmin=None, bc_zmax=None, moving_window_velocity=None, refined_regions=[], lower_bound_particles=None, upper_bound_particles=None, rmin_particles=None, rmax_particles=None, zmin_particles=None, zmax_particles=None, lower_boundary_conditions_particles=None, upper_boundary_conditions_particles=None, bc_rmin_particles=None, bc_rmax_particles=None, bc_zmin_particles=None, bc_zmax_particles=None, guard_cells=None, pml_cells=None, **kw)[source]

Axisymmetric, cylindrical grid Parameters can be specified either as vectors or separately. (If both are specified, the vector is used.)

Parameters
  • number_of_cells (vector of integers) – Number of cells along each axis (number of nodes is number_of_cells+1)

  • lower_bound (vector of floats) – Position of the node at the lower bound [m]

  • upper_bound (vector of floats) – Position of the node at the upper bound [m]

  • lower_boundary_conditions (vector of strings) – Conditions at lower boundaries, periodic, open, dirichlet, or neumann

  • upper_boundary_conditions (vector of strings) – Conditions at upper boundaries, periodic, open, dirichlet, or neumann

  • nr (integer) – Number of cells along R (number of nodes=nr+1)

  • nz (integer) – Number of cells along Z (number of nodes=nz+1)

  • n_azimuthal_modes (integer) – Number of azimuthal modes

  • rmin (float) – Position of first node along R [m]

  • rmax (float) – Position of last node along R [m]

  • zmin (float) – Position of first node along Z [m]

  • zmax (float) – Position of last node along Z [m]

  • bc_rmin (vector of strings) – Boundary condition at min R: One of open, dirichlet, or neumann

  • bc_rmax (vector of strings) – Boundary condition at max R: One of open, dirichlet, or neumann

  • bc_zmin (vector of strings) – Boundary condition at min Z: One of periodic, open, dirichlet, or neumann

  • bc_zmax (vector of strings) – Boundary condition at max Z: One of periodic, open, dirichlet, or neumann

  • moving_window_velocity (vector of floats, optional) – Moving frame velocity [m/s]

  • refined_regions (list of lists, optional) – List of refined regions, each element being a list of the format [level, lo, hi, refinement_factor], with level being the refinement level, with 1 being the first level of refinement, 2 being the second etc, lo and hi being vectors of length 2 specifying the extent of the region, and refinement_factor defaulting to [2,2] (relative to next lower level)

  • lower_bound_particles (vector of floats, optional) – Position of particle lower bound [m]

  • upper_bound_particles (vector of floats, optional) – Position of particle upper bound [m]

  • rmin_particles (float, optional) – Position of min particle boundary along R [m]

  • rmax_particles (float, optional) – Position of max particle boundary along R [m]

  • zmin_particles (float, optional) – Position of min particle boundary along Z [m]

  • zmax_particles (float, optional) – Position of max particle boundary along Z [m]

  • lower_boundary_conditions_particles (vector of strings, optional) – Conditions at lower boundaries for particles, periodic, absorbing, reflect or thermal

  • upper_boundary_conditions_particles (vector of strings, optional) – Conditions at upper boundaries for particles, periodic, absorbing, reflect or thermal

  • bc_rmin_particles (string, optional) – Boundary condition at min R for particles: One of periodic, absorbing, reflect, thermal

  • bc_rmax_particles (string, optional) – Boundary condition at max R for particles: One of periodic, absorbing, reflect, thermal

  • bc_zmin_particles (string, optional) – Boundary condition at min Z for particles: One of periodic, absorbing, reflect, thermal

  • bc_zmax_particles (string, optional) – Boundary condition at max Z for particles: One of periodic, absorbing, reflect, thermal

  • guard_cells (vector of integers, optional) – Number of guard cells used along each direction

  • pml_cells (vector of integers, optional) – Number of Perfectly Matched Layer (PML) cells along each direction

Implementation specific documentation

This assumes that WarpX was compiled with USE_RZ = TRUE

See Input Parameters for more information.

Parameters
  • warpx_max_grid_size (integer, default=32) – Maximum block size in either direction

  • warpx_max_grid_size_x (integer, optional) – Maximum block size in radial direction

  • warpx_max_grid_size_y (integer, optional) – Maximum block size in longitudinal direction

  • warpx_blocking_factor (integer, optional) – Blocking factor (which controls the block size)

  • warpx_blocking_factor_x (integer, optional) – Blocking factor (which controls the block size) in the radial direction

  • warpx_blocking_factor_y (integer, optional) – Blocking factor (which controls the block size) in the longitudinal direction

  • warpx_potential_lo_r (float, default=0.) – Electrostatic potential on the lower radial boundary

  • warpx_potential_hi_r (float, default=0.) – Electrostatic potential on the upper radial boundary

  • warpx_potential_lo_z (float, default=0.) – Electrostatic potential on the lower longitudinal boundary

  • warpx_potential_hi_z (float, default=0.) – Electrostatic potential on the upper longitudinal boundary

EmbeddedBoundary

class pywarpx.picmi.EmbeddedBoundary(implicit_function=None, stl_file=None, stl_scale=None, stl_center=None, stl_reverse_normal=False, potential=None, **kw)[source]

Custom class to handle set up of embedded boundaries specific to WarpX. If embedded boundary initialization is added to picmistandard this can be changed to inherit that functionality. The geometry can be specified either as an implicit function or as an STL file (ASCII or binary). In the latter case the geometry specified in the STL file can be scaled, translated and inverted.

Parameters
  • implicit_function (string) – Analytic expression describing the embedded boundary

  • stl_file (string) – STL file path (string), file contains the embedded boundary geometry

  • stl_scale (float) – Factor by which the STL geometry is scaled

  • stl_center (vector of floats) – Vector by which the STL geometry is translated (in meters)

  • stl_reverse_normal (bool) – If True inverts the orientation of the STL geometry

  • potential (string, default=0.) – Analytic expression defining the potential. Can only be specified when the solver is electrostatic.

Parameters used in the analytic expressions should be given as additional keyword arguments.

Applied fields

ConstantAppliedField

class pywarpx.picmi.ConstantAppliedField(Ex=None, Ey=None, Ez=None, Bx=None, By=None, Bz=None, lower_bound=[None, None, None], upper_bound=[None, None, None], **kw)[source]

Describes a constant applied field

Parameters
  • Ex (float, default=0.) – Constant Ex field [V/m]

  • Ey (float, default=0.) – Constant Ey field [V/m]

  • Ez (float, default=0.) – Constant Ez field [V/m]

  • Bx (float, default=0.) – Constant Bx field [T]

  • By (float, default=0.) – Constant By field [T]

  • Bz (float, default=0.) – Constant Bz field [T]

  • lower_bound (vector, optional) – Lower bound of the region where the field is applied [m].

  • upper_bound (vector, optional) – Upper bound of the region where the field is applied [m]

AnalyticAppliedField

class pywarpx.picmi.AnalyticAppliedField(Ex_expression=None, Ey_expression=None, Ez_expression=None, Bx_expression=None, By_expression=None, Bz_expression=None, lower_bound=[None, None, None], upper_bound=[None, None, None], **kw)[source]

Describes an analytic applied field

The expressions should be in terms of the position and time, written as ‘x’, ‘y’, ‘z’, ‘t’. Parameters can be used in the expression with the values given as additional keyword arguments. Expressions should be relative to the lab frame.

Parameters
  • Ex_expression (string, optional) – Analytic expression describing Ex field [V/m]

  • Ey_expression (string, optional) – Analytic expression describing Ey field [V/m]

  • Ez_expression (string, optional) – Analytic expression describing Ez field [V/m]

  • Bx_expression (string, optional) – Analytic expression describing Bx field [T]

  • By_expression (string, optional) – Analytic expression describing By field [T]

  • Bz_expression (string, optional) – Analytic expression describing Bz field [T]

  • lower_bound (vector, optional) – Lower bound of the region where the field is applied [m].

  • upper_bound (vector, optional) – Upper bound of the region where the field is applied [m]

PlasmaLens

class pywarpx.picmi.PlasmaLens(period, starts, lengths, strengths_E=None, strengths_B=None, **kw)[source]

Custom class to setup a plasma lens lattice. The applied fields are dependent only on the transverse position.

Parameters
  • period (float) – Periodicity of the lattice (in lab frame, in meters)

  • starts (list of floats) – The start of each lens relative to the periodic repeat

  • lengths (list of floats) – The length of each lens

  • strengths_E=None (list of floats, default = 0.) – The electric field strength of each lens

  • strengths_B=None (list of floats, default = 0.) – The magnetic field strength of each lens

The field that is applied depends on the transverse position of the particle, (x,y)

  • Ex = x*stengths_E

  • Ey = y*stengths_E

  • Bx = +y*stengths_B

  • By = -x*stengths_B

Mirror

class pywarpx.picmi.Mirror(x_front_location=None, y_front_location=None, z_front_location=None, depth=None, number_of_cells=None, **kw)[source]

Describes a perfectly reflecting mirror, where the E and B fields are zeroed out in a plane of finite thickness.

Parameters
  • x_front_location (float, optional (see comment below)) – Location in x of the front of the nirror [m]

  • y_front_location (float, optional (see comment below)) – Location in y of the front of the nirror [m]

  • z_front_location (float, optional (see comment below)) – Location in z of the front of the nirror [m]

  • depth (float, optional (see comment below)) – Depth of the mirror [m]

  • number_of_cells (integer, optional (see comment below)) – Minimum numer of cells zeroed out

Only one of the [x,y,z]_front_location should be specified. The mirror will be set perpendicular to the respective direction and infinite in the others. The depth of the mirror will be the maximum of the specified depth and number_of_cells, or the code’s default value if neither are specified.

Diagnostics

ParticleDiagnostic

class pywarpx.picmi.ParticleDiagnostic(period, species, data_list=None, write_dir=None, step_min=None, step_max=None, parallelio=None, name=None, **kw)[source]

Defines the particle diagnostics in the simulation frame

Parameters
  • period (integer) – Period of time steps that the diagnostic is performed

  • species (species instance or list of species instances) – Species to write out. Note that the name attribute must be defined for the species.

  • data_list (list of strings, optional) – The data to be written out. Possible values ‘position’, ‘momentum’, ‘weighting’. Defaults to the output list of the implementing code.

  • write_dir (string, optional) – Directory where data is to be written

  • step_min (integer, default=0) – Minimum step at which diagnostics could be written

  • step_max (integer, default=unbounded) – Maximum step at which diagnostics could be written

  • parallelio (bool, optional) – If set to True, particle diagnostics are dumped in parallel

  • name (string, optional) – Sets the base name for the diagnostic output files

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_format ({plotfile, checkpoint, openpmd, ascent, sensei}, optional) – Diagnostic file format

  • warpx_openpmd_backend ({bp, h5, json}, optional) – Openpmd backend file format

  • warpx_file_prefix (string, optional) – Prefix on the diagnostic file name

  • warpx_file_min_digits (integer, optional) – Minimum number of digits for the time step number in the file name

  • warpx_random_fraction (float, optional) – Random fraction of particles to include in the diagnostic

  • warpx_uniform_stride (integer, optional) – Stride to down select to the particles to include in the diagnostic

  • warpx_plot_filter_function (string, optional) – Analytic expression to down select the particles to in the diagnostic

FieldDiagnostic

class pywarpx.picmi.FieldDiagnostic(grid, period, data_list=None, write_dir=None, step_min=None, step_max=None, number_of_cells=None, lower_bound=None, upper_bound=None, parallelio=None, name=None, **kw)[source]

Defines the electromagnetic field diagnostics in the simulation frame

Parameters
  • grid (grid instance) – Grid object for the diagnostic

  • period (integer) – Period of time steps that the diagnostic is performed

  • data_list (list of strings, optional) – List of quantities to write out. Possible values ‘rho’, ‘E’, ‘B’, ‘J’, ‘Ex’ etc. Defaults to the output list of the implementing code.

  • write_dir (string, optional) – Directory where data is to be written

  • step_min (integer, default=0) – Minimum step at which diagnostics could be written

  • step_max (integer, default=unbounded) – Maximum step at which diagnostics could be written

  • number_of_cells (vector of integers, optional) – Number of cells in each dimension. If not given, will be obtained from grid.

  • lower_bound (vector of floats, optional) – Lower corner of diagnostics box in each direction. If not given, will be obtained from grid.

  • upper_bound (vector of floats, optional) – Higher corner of diagnostics box in each direction. If not given, will be obtained from grid.

  • parallelio (bool, optional) – If set to True, field diagnostics are dumped in parallel

  • name (string, optional) – Sets the base name for the diagnostic output files

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_plot_raw_fields (bool, optional) – Flag whether to dump the raw fields

  • warpx_plot_raw_fields_guards (bool, optional) – Flag whether the raw fields should include the guard cells

  • warpx_format ({plotfile, checkpoint, openpmd, ascent, sensei}, optional) – Diagnostic file format

  • warpx_openpmd_backend ({bp, h5, json}, optional) – Openpmd backend file format

  • warpx_file_prefix (string, optional) – Prefix on the diagnostic file name

  • warpx_file_min_digits (integer, optional) – Minimum number of digits for the time step number in the file name

  • warpx_dump_rz_modes (bool, optional) – Flag whether to dump the data for all RZ modes

ElectrostaticFieldDiagnostic

pywarpx.picmi.ElectrostaticFieldDiagnostic

alias of FieldDiagnostic

Lab-frame diagnostics diagnostics are used when running boosted-frame simulations.

LabFrameParticleDiagnostic

class pywarpx.picmi.LabFrameParticleDiagnostic(grid, num_snapshots, dt_snapshots, data_list=None, time_start=0.0, species=None, write_dir=None, parallelio=None, name=None, **kw)[source]

Defines the particle diagnostics in the lab frame

Parameters
  • grid (grid instance) – Grid object for the diagnostic

  • num_snapshots (integer) – Number of lab frame snapshots to make

  • dt_snapshots (float) – Time between each snapshot in lab frame

  • species (species instance or list of species instances) – Species to write out. Note that the name attribute must be defined for the species.

  • data_list (list of strings, optional) – The data to be written out. Possible values ‘position’, ‘momentum’, ‘weighting’. Defaults to the output list of the implementing code.

  • time_start (float, default=0) – Time for the first snapshot in lab frame

  • write_dir (string, optional) – Directory where data is to be written

  • parallelio (bool, optional) – If set to True, particle diagnostics are dumped in parallel

  • name (string, optional) – Sets the base name for the diagnostic output files

LabFrameFieldDiagnostic

class pywarpx.picmi.LabFrameFieldDiagnostic(grid, num_snapshots, dt_snapshots, data_list=None, z_subsampling=1, time_start=0.0, write_dir=None, parallelio=None, name=None, **kw)[source]

Defines the electromagnetic field diagnostics in the lab frame

Parameters
  • grid (grid instance) – Grid object for the diagnostic

  • num_snapshots (integer) – Number of lab frame snapshots to make

  • dt_snapshots (float) – Time between each snapshot in lab frame

  • data_list (list of strings, optional) – List of quantities to write out. Possible values ‘rho’, ‘E’, ‘B’, ‘J’, ‘Ex’ etc. Defaults to the output list of the implementing code.

  • z_subsampling (integer, default=1) – A factor which is applied on the resolution of the lab frame reconstruction

  • time_start (float, default=0) – Time for the first snapshot in lab frame

  • write_dir (string, optional) – Directory where data is to be written

  • parallelio (bool, optional) – If set to True, field diagnostics are dumped in parallel

  • name (string, optional) – Sets the base name for the diagnostic output files

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_new_BTD (bool, optional) – Use the new BTD diagnostics

  • warpx_format (string, optional) – Passed to <diagnostic name>.format

  • warpx_openpmd_backend (string, optional) – Passed to <diagnostic name>.openpmd_backend

  • warpx_file_prefix (string, optional) – Passed to <diagnostic name>.file_prefix

  • warpx_file_min_digits (integer, optional) – Passed to <diagnostic name>.file_min_digits

  • warpx_buffer_size (integer, optional) – Passed to <diagnostic name>.buffer_size

  • warpx_lower_bound (vector of floats, optional) – Passed to <diagnostic name>.lower_bound

  • warpx_upper_bound (vector of floats, optional) – Passed to <diagnostic name>.upper_bound

Checkpoint

class pywarpx.picmi.Checkpoint(period=1, write_dir=None, name=None, **kw)[source]

Sets up checkpointing of the simulation, allowing for later restarts

See Input Parameters for more information.

Parameters
  • warpx_file_prefix (string) – The prefix to the checkpoint directory names

  • warpx_file_min_digits (integer) – Minimum number of digits for the time step number in the checkpoint directory name.

Particles

Species objects are a collection of particles with similar properties. For instance, background plasma electrons, background plasma ions and an externally injected beam could each be their own particle species.

Species

class pywarpx.picmi.Species(particle_type=None, name=None, charge_state=None, charge=None, mass=None, initial_distribution=None, particle_shape=None, density_scale=None, method=None, **kw)[source]

Sets up the species to be simulated. The species charge and mass can be specified by setting the particle type or by setting them directly. If the particle type is specified, the charge or mass can be set to override the value from the type.

Parameters
  • particle_type (string, optional) – A string specifying an elementary particle, atom, or other, as defined in the openPMD 2 species type extension, openPMD-standard/EXT_SpeciesType.md

  • name (string, optional) – Name of the species

  • method ({'Boris', 'Vay', 'Higuera-Cary', 'Li' , 'free-streaming', 'LLRK4'}) –

    The particle advance method to use. Code-specific method can be specified using ‘other:<method>’. The default is code dependent.

    • ’Boris’: Standard “leap-frog” Boris advance

    • ’Vay’:

    • ’Higuera-Cary’:

    • ’Li’ :

    • ’free-streaming’: Advance with no fields

    • ’LLRK4’: Landau-Lifschitz radiation reaction formula with RK-4)

  • charge_state (float, optional) – Charge state of the species (applies only to atoms) [1]

  • charge (float, optional) – Particle charge, required when type is not specified, otherwise determined from type [C]

  • mass (float, optional) – Particle mass, required when type is not specified, otherwise determined from type [kg]

  • initial_distribution (distribution instance) – The initial distribution loaded at t=0. Must be one of the standard distributions implemented.

  • density_scale (float, optional) – A scale factor on the density given by the initial_distribution

  • particle_shape ({'NGP', 'linear', 'quadratic', 'cubic'}) – Particle shape used for deposition and gather. If not specified, the value from the Simulation object will be used. Other values maybe specified that are code dependent.

Implementation specific documentation

See Input Parameters for more information.

Parameters
  • warpx_boost_adjust_transverse_positions (bool, default=False) – Whether to adjust transverse positions when apply the boost to the simulation frame

  • warpx_self_fields_required_precision (float, default=1.e-11) – Relative precision on the electrostatic solver (when using the relativistic solver)

  • warpx_self_fields_absolute_tolerance (float, default=0.) – Absolute precision on the electrostatic solver (when using the relativistic solver)

  • warpx_self_fields_max_iters (integer, default=200) – Maximum number of iterations for the electrostatic solver for the species

  • warpx_self_fields_verbosity (integer, default=2) – Level of verbosity for the electrostatic solver

  • warpx_save_previous_position (bool, default=False) – Whether to save the old particle positions

  • warpx_do_not_deposit (bool, default=False) – Whether or not to deposit the charge and current density for for this species

  • warpx_reflection_model_xlo (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the lower x boundary

  • warpx_reflection_model_xhi (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the upper x boundary

  • warpx_reflection_model_ylo (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the lower y boundary

  • warpx_reflection_model_yhi (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the upper y boundary

  • warpx_reflection_model_zlo (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the lower z boundary

  • warpx_reflection_model_zhi (string, default='0.') – Expression (in terms of the velocity “v”) specifying the probability that the particle will reflect on the upper z boundary

  • warpx_save_particles_at_xlo (bool, default=False) – Whether to save particles lost at the lower x boundary

  • warpx_save_particles_at_xhi (bool, default=False) – Whether to save particles lost at the upper x boundary

  • warpx_save_particles_at_ylo (bool, default=False) – Whether to save particles lost at the lower y boundary

  • warpx_save_particles_at_yhi (bool, default=False) – Whether to save particles lost at the upper y boundary

  • warpx_save_particles_at_zlo (bool, default=False) – Whether to save particles lost at the lower z boundary

  • warpx_save_particles_at_zhi (bool, default=False) – Whether to save particles lost at the upper z boundary

  • warpx_save_particles_at_eb (bool, default=False) – Whether to save particles lost at the embedded boundary

MultiSpecies

class pywarpx.picmi.MultiSpecies(particle_types=None, names=None, charge_states=None, charges=None, masses=None, proportions=None, initial_distribution=None, particle_shape=None, **kw)[source]

INCOMPLETE: proportions argument is not implemented Multiple species that are initialized with the same distribution. Each parameter can be list, giving a value for each species, or a single value which is given to all species. The species charge and mass can be specified by setting the particle type or by setting them directly. If the particle type is specified, the charge or mass can be set to override the value from the type.

Parameters
  • particle_types (list of strings, optional) – A string specifying an elementary particle, atom, or other, as defined in the openPMD 2 species type extension, openPMD-standard/EXT_SpeciesType.md

  • names (list of strings, optional) – Names of the species

  • charge_states (list of floats, optional) – Charge states of the species (applies only to atoms)

  • charges (list of floats, optional) – Particle charges, required when type is not specified, otherwise determined from type [C]

  • masses (list of floats, optional) – Particle masses, required when type is not specified, otherwise determined from type [kg]

  • proportions (list of floats, optional) – Proportions of the initial distribution made up by each species

  • initial_distribution (distribution instance) – Initial particle distribution, applied to all species

  • particle_shape ({'NGP', 'linear', 'quadratic', 'cubic'}) – Particle shape used for deposition and gather. If not specified, the value from the Simulation object will be used. Other values maybe specified that are code dependent.

Particle distributions can be used for to initialize particles in a particle species.

GaussianBunchDistribution

class pywarpx.picmi.GaussianBunchDistribution(n_physical_particles, rms_bunch_size, rms_velocity=[0.0, 0.0, 0.0], centroid_position=[0.0, 0.0, 0.0], centroid_velocity=[0.0, 0.0, 0.0], velocity_divergence=[0.0, 0.0, 0.0], **kw)[source]

Describes a Gaussian distribution of particles

Parameters
  • n_physical_particles (integer) – Number of physical particles in the bunch

  • rms_bunch_size (vector of floats) – RMS bunch size at t=0 [m]

  • rms_velocity (vector of floats, default=[0.,0.,0.]) – RMS velocity spread at t=0 [m/s]

  • centroid_position (vector of floats, default=[0.,0.,0.]) – Position of the bunch centroid at t=0 [m]

  • centroid_velocity (vector of floats, default=[0.,0.,0.]) – Velocity (gamma*V) of the bunch centroid at t=0 [m/s]

  • velocity_divergence (vector of floats, default=[0.,0.,0.]) – Expansion rate of the bunch at t=0 [m/s/m]

UniformDistribution

class pywarpx.picmi.UniformDistribution(density, lower_bound=[None, None, None], upper_bound=[None, None, None], rms_velocity=[0.0, 0.0, 0.0], directed_velocity=[0.0, 0.0, 0.0], fill_in=None, **kw)[source]

Describes a uniform density distribution of particles

Parameters
  • density (float) – Physical number density [m^-3]

  • lower_bound (vector of floats, optional) – Lower bound of the distribution [m]

  • upper_bound (vector of floats, optional) – Upper bound of the distribution [m]

  • rms_velocity (vector of floats, default=[0.,0.,0.]) – Thermal velocity spread [m/s]

  • directed_velocity (vector of floats, default=[0.,0.,0.]) – Directed, average, velocity [m/s]

  • fill_in (bool, optional) – Flags whether to fill in the empty spaced opened up when the grid moves

AnalyticDistribution

class pywarpx.picmi.AnalyticDistribution(density_expression, momentum_expressions=[None, None, None], lower_bound=[None, None, None], upper_bound=[None, None, None], rms_velocity=[0.0, 0.0, 0.0], directed_velocity=[0.0, 0.0, 0.0], fill_in=None, **kw)[source]

Describes a uniform density plasma

Parameters
  • density_expression (string) – Analytic expression describing physical number density (string) [m^-3]. Expression should be in terms of the position, written as ‘x’, ‘y’, and ‘z’. Parameters can be used in the expression with the values given as keyword arguments.

  • momentum_expressions (list of strings) – Analytic expressions describing the gamma*velocity for each axis [m/s]. Expressions should be in terms of the position, written as ‘x’, ‘y’, and ‘z’. Parameters can be used in the expression with the values given as keyword arguments. For any axis not supplied (set to None), directed_velocity will be used.

  • lower_bound (vector of floats, optional) – Lower bound of the distribution [m]

  • upper_bound (vector of floats, optional) – Upper bound of the distribution [m]

  • rms_velocity (vector of floats, detault=[0.,0.,0.]) – Thermal velocity spread [m/s]

  • directed_velocity (vector of floats, detault=[0.,0.,0.]) – Directed, average, velocity [m/s]

  • fill_in (bool, optional) – Flags whether to fill in the empty spaced opened up when the grid moves

This will create a distribution where the density is n0 below rmax and zero elsewhere.:

.. code-block: python
dist = AnalyticDistribution(density_expression=’((x**2+y**2)<rmax**2)*n0’,

rmax = 1., n0 = 1.e20, …)

ParticleListDistribution

class pywarpx.picmi.ParticleListDistribution(x=0.0, y=0.0, z=0.0, ux=0.0, uy=0.0, uz=0.0, weight=0.0, **kw)[source]

Load particles at the specified positions and velocities

Parameters
  • x (float, default=0.) – List of x positions of the particles [m]

  • y (float, default=0.) – List of y positions of the particles [m]

  • z (float, default=0.) – List of z positions of the particles [m]

  • ux (float, default=0.) – List of ux positions of the particles (ux = gamma*vx) [m/s]

  • uy (float, default=0.) – List of uy positions of the particles (uy = gamma*vy) [m/s]

  • uz (float, default=0.) – List of uz positions of the particles (uz = gamma*vz) [m/s]

  • weight (float) – Particle weight or list of weights, number of real particles per simulation particle

Particle layouts determine how to microscopically place macro particles in a grid cell.

GriddedLayout

class pywarpx.picmi.GriddedLayout(n_macroparticle_per_cell, grid=None, **kw)[source]

Specifies a gridded layout of particles

Parameters
  • n_macroparticle_per_cell (vector of integers) – Number of particles per cell along each axis

  • grid (grid instance, optional) – Grid object specifying the grid to follow. If not specified, the underlying grid of the code is used.

PseudoRandomLayout

class pywarpx.picmi.PseudoRandomLayout(n_macroparticles=None, n_macroparticles_per_cell=None, seed=None, grid=None, **kw)[source]

Specifies a pseudo-random layout of the particles

Parameters
  • n_macroparticles (integer) – Total number of macroparticles to load. Either this argument or n_macroparticles_per_cell should be supplied.

  • n_macroparticles_per_cell (integer) – Number of macroparticles to load per cell. Either this argument or n_macroparticles should be supplied.

  • seed (integer, optional) – Pseudo-random number generator seed

  • grid (grid instance, optional) – Grid object specifying the grid to follow for n_macroparticles_per_cell. If not specified, the underlying grid of the code is used.

Other operations related to particles

CoulombCollisions

class pywarpx.picmi.CoulombCollisions(name, species, CoulombLog=None, ndt=None, **kw)[source]

Custom class to handle setup of binary Coulmb collisions in WarpX. If collision initialization is added to picmistandard this can be changed to inherit that functionality.

Parameters
  • name (string) – Name of instance (used in the inputs file)

  • species (list of species instances) – The species involved in the collision. Must be of length 2.

  • CoulombLog (float, optional) – Value of the Coulomb log to use in the collision cross section. If not supplied, it is calculated from the local conditions.

  • ndt (integer, optional) – The collisions will be applied every “ndt” steps. Must be 1 or larger.

MCCCollisions

class pywarpx.picmi.MCCCollisions(name, species, background_density, background_temperature, scattering_processes, background_mass=None, ndt=None, **kw)[source]

Custom class to handle setup of MCC collisions in WarpX. If collision initialization is added to picmistandard this can be changed to inherit that functionality.

Parameters
  • name (string) – Name of instance (used in the inputs file)

  • species (species instance) – The species involved in the collision

  • background_density (float) – The density of the background

  • background_temperature (float) – The temperature of the background

  • scattering_processes (dictionary) – The scattering process to use and any needed information

  • background_mass (float, optional) – The mass of the background particle. If not supplied, the default depends on the type of scattering process.

  • ndt (integer, optional) – The collisions will be applied every “ndt” steps. Must be 1 or larger.

Lasers

Laser profiles can be used to initialize laser pulses in the simulation.

GaussianLaser

class pywarpx.picmi.GaussianLaser(wavelength, waist, duration, focal_position=[0.0, 0.0, 0.0], centroid_position=[0.0, 0.0, 0.0], propagation_direction=[0.0, 0.0, 1.0], polarization_direction=[1.0, 0.0, 0.0], a0=None, E0=None, phi0=None, zeta=None, beta=None, phi2=None, name=None, fill_in=True, **kw)[source]

Specifies a Gaussian laser distribution.

More precisely, the electric field near the focal plane is given by:

\[E(\boldsymbol{x},t) = a_0\times E_0\, \exp\left( -\frac{r^2}{w_0^2} - \frac{(z-z_0-ct)^2}{c^2\tau^2} \right) \cos[ k_0( z - z_0 - ct ) - \phi_{cep} ]\]

where \(k_0 = 2\pi/\lambda_0\) is the wavevector and where \(E_0 = m_e c^2 k_0 / q_e\) is the field amplitude for \(a_0=1\).

Note

The additional terms that arise far from the focal plane (Gouy phase, wavefront curvature, …) are not included in the above formula for simplicity, but are of course taken into account by the code, when initializing the laser pulse away from the focal plane.

Parameters
  • wavelength (float) – Laser wavelength [m], defined as \(\lambda_0\) in the above formula

  • waist (float) – Waist of the Gaussian pulse at focus [m], defined as \(w_0\) in the above formula

  • duration (float) – Duration of the Gaussian pulse [s], defined as \(\tau\) in the above formula

  • focal_position (vector of floats, default=[0,0,0]) – Position of the laser focus [m]

  • centroid_position (vector of floats, default=[0,0,0]) – Position of the laser centroid at time 0 [m]

  • propagation_direction (unit vector of floats, default=[0,0,1]) – Direction of propagation [1]

  • polarization_direction (unit vector of floats, default=[1,0,0]) – Direction of polarization [1]

  • a0 (float) – Normalized vector potential at focus Specify either a0 or E0 (E0 takes precedence).

  • E0 (float) – Maximum amplitude of the laser field [V/m] Specify either a0 or E0 (E0 takes precedence).

  • phi0 (float) – Carrier envelope phase (CEP) [rad]

  • zeta (float) – Spatial chirp at focus (in the lab frame) [m.s]

  • beta (float) – Angular dispersion at focus (in the lab frame) [rad.s]

  • phi2 (float) – Temporal chirp at focus (in the lab frame) [s^2]

  • fill_in (bool, default=True) – Flags whether to fill in the empty spaced opened up when the grid moves

  • name (string, optional) – Optional name of the laser

AnalyticLaser

class pywarpx.picmi.AnalyticLaser(field_expression, wavelength, propagation_direction=[0.0, 0.0, 1.0], polarization_direction=[1.0, 0.0, 0.0], amax=None, Emax=None, name=None, fill_in=True, **kw)[source]

Specifies a laser with an analytically described distribution

Parameters
  • name=None (string, optional) – Optional name of the laser

  • field_expression (string) – Analytic expression describing the electric field of the laser [V/m] Expression should be in terms of the position, ‘X’, ‘Y’, in the plane orthogonal to the propagation direction, and ‘t’ the time. The expression should describe the full field, including the oscillitory component. Parameters can be used in the expression with the values given as keyword arguments.

  • propagation_direction (unit vector of floats, default=[0,0,1]) – Direction of propagation [1]

  • polarization_direction (unit vector of floats, default=[1,0,0]) – Direction of polarization [1]

  • wavelength (float, optional) – Laser wavelength. This should be built into the expression, but some codes require a specified value for numerical purposes.

  • amax (float, optional) – Maximum normalized vector potential. Specify either amax or Emax (Emax takes precedence). This should be built into the expression, but some codes require a specified value for numerical purposes.

  • Emax (float, optional) – Maximum amplitude of the laser field [V/m]. Specify either amax or Emax (Emax takes precedence). This should be built into the expression, but some codes require a specified value for numerical purposes.

  • fill_in (bool, default=True) – Flags whether to fill in the empty spaced opened up when the grid moves

Laser injectors control where to initialize laser pulses on the simulation grid.

LaserAntenna

class pywarpx.picmi.LaserAntenna(position, normal_vector, **kw)[source]

Specifies the laser antenna injection method

Parameters
  • position (vector of strings) – Position of antenna launching the laser [m]

  • normal_vector (vector of strings) – Vector normal to antenna plane [1]

Running

WarpX can be run in one of two modes. It can run as a preprocessor, using the Python input file to generate an input file to be used by the C++ version, or it can be run directly from Python. The examples support running in both modes by commenting and uncommenting the appropriate lines.

In either mode, if using a virtual environment, be sure to activate it before compiling and running WarpX.

Running WarpX directly from Python

For this, a full Python installation of WarpX is required, as described in the install documentation (developers).

In order to run a new simulation:

  • Create a new directory, where the simulation will be run.

  • Add a Python script in the directory.

The input file should have the line sim.step() which runs the simulation.

  • Run the script with Python:

mpirun -np <n_ranks> python <python_script>

where <n_ranks> is the number of MPI ranks used, and <python_script> is the name of the script.

Extending a Simulation from Python

When running WarpX directly from Python it is possible to interact with the simulation by installing CallbackFunctions, which will execute a given Python function at a specific location in the WarpX simulation loop.

class pywarpx.callbacks.CallbackFunctions(name=None, lcallonce=0)[source]

Class to handle call back function lists.

Note that for functions passed in that are methods of a class instance, a full reference of the instance is saved. This extra reference means that the object will not actually be deleted if the user deletes the original reference. This is good since the user does not need to keep the reference to the object (for example it can be created using a local variable in a function). It may be bad if the user thinks an object was deleted, but it actually isn’t since it had (unkown to the user) installed a method in one of the call back lists.

Places in the WarpX loop where callbacks are available include: afterinit, beforecollisions, aftercollisions, beforeEsolve, afterEsolve, beforedeposition, afterdeposition, beforestep, afterstep, afterdiagnostics, afterrestart and oncheckpointsignal. See the examples in Examples/Tests/ParticleDataPython for references on how to use callbacks.

There are several “hooks” available via the libwarpx shared library to access and manipulate simulation objects (particles, fields and memory buffers) as well as general properties (such as processor number). These “hooks” are accessible through the Simulation.extension object.

pywarpx.picmi.Simulation.extension.getNProcs()

Get the number of processors

pywarpx.picmi.Simulation.extension.getMyProc()

Get the number of the processor

pywarpx.picmi.Simulation.extension.get_nattr()

Get the number of extra particle attributes.

pywarpx.picmi.Simulation.extension.get_nattr_species(species_name)

Get the number of real attributes for the given species.

Parameters

species_name (str) – Name of the species

pywarpx.picmi.Simulation.extension.getistep(level=0)

Get the current time step number for the specified level

Parameter

levelint

The refinement level to reference

pywarpx.picmi.Simulation.extension.gett_new(level=0)

Get the next time for the specified level.

Parameters

level (int) – The refinement level to reference

pywarpx.picmi.Simulation.extension.evolve(num_steps=-1)

Evolve the simulation for num_steps steps. If num_steps=-1, the simulation will be run until the end as specified in the inputs file.

Parameters

num_steps (int) – The number of steps to take

pywarpx.picmi.Simulation.extension.finalize(finalize_mpi=1)

Call finalize for WarpX and AMReX. Registered to run at program exit.

pywarpx.picmi.Simulation.extension.getistep(level=0)

Get the current time step number for the specified level

Parameter

levelint

The refinement level to reference

pywarpx.picmi.Simulation.extension.gett_new(level=0)

Get the next time for the specified level.

Parameters

level (int) – The refinement level to reference

pywarpx.picmi.Simulation.extension.evolve(num_steps=-1)

Evolve the simulation for num_steps steps. If num_steps=-1, the simulation will be run until the end as specified in the inputs file.

Parameters

num_steps (int) – The number of steps to take

pywarpx.picmi.Simulation.extension.getProbLo(direction)

Get the values of the lower domain boundary.

Parameters

direction (int) – Direction of interest

pywarpx.picmi.Simulation.extension.getProbHi(direction)

Get the values of the upper domain boundary.

Parameters

direction (int) – Direction of interest

pywarpx.picmi.Simulation.extension.getCellSize(direction, level=0)

Get the cell size in the given direction and on the given level.

Parameters
  • direction (int) – Direction of interest

  • level (int) – The refinement level to reference

Particles can be added to the simulation at specific positions and with specific attribute values:

pywarpx.picmi.Simulation.extension.add_particles(species_name, x=None, y=None, z=None, ux=None, uy=None, uz=None, w=None, unique_particles=True, **kwargs)

A function for adding particles to the WarpX simulation.

Parameters
  • species_name (str) – The type of species for which particles will be added

  • x (arrays or scalars) – The particle positions (default = 0.)

  • y (arrays or scalars) – The particle positions (default = 0.)

  • z (arrays or scalars) – The particle positions (default = 0.)

  • ux (arrays or scalars) – The particle momenta (default = 0.)

  • uy (arrays or scalars) – The particle momenta (default = 0.)

  • uz (arrays or scalars) – The particle momenta (default = 0.)

  • w (array or scalars) – Particle weights (default = 0.)

  • unique_particles (bool) – Whether the particles are unique or duplicated on several processes (default = True)

  • kwargs (dict) – Containing an entry for all the extra particle attribute arrays. If an attribute is not given it will be set to 0.

Properties of the particles already in the simulation can be obtained with various functions.

pywarpx.picmi.Simulation.extension.get_particle_count(species_name, local=False)

Get the number of particles of the specified species in the simulation.

Parameters
  • species_name (str) – The species name that the number will be returned for

  • local (bool) – If True the particle count on this processor will be returned. Default False.

Returns

An integer count of the number of particles

Return type

int

pywarpx.picmi.Simulation.extension.get_particle_structs(species_name, level)

This returns a list of numpy arrays containing the particle struct data on each tile for this process. The particle data is represented as a structured numpy array and contains the particle ‘x’, ‘y’, ‘z’, and ‘idcpu’.

The data for the numpy arrays are not copied, but share the underlying memory buffer with WarpX. The numpy arrays are fully writeable.

Parameters
  • species_name (str) – The species name that the data will be returned for

  • level (int) – The refinement level to reference

Returns

The requested particle struct data

Return type

List of numpy arrays

pywarpx.picmi.Simulation.extension.get_particle_arrays(species_name, comp_name, level)

This returns a list of numpy arrays containing the particle array data on each tile for this process.

The data for the numpy arrays are not copied, but share the underlying memory buffer with WarpX. The numpy arrays are fully writeable.

Parameters
  • species_name (str) – The species name that the data will be returned for

  • comp_name (str) – The component of the array data that will be returned

  • level (int) – The refinement level to reference

Returns

The requested particle array data

Return type

List of numpy arrays

The get_particle_structs() and get_particle_arrays() functions are called by several utility functions of the form get_particle_{comp_name} where comp_name is one of x, y, z, r, theta, id, cpu, weight, ux, uy or uz.

The index of some specific component of the particle data can be obtained.

pywarpx.picmi.Simulation.extension.get_particle_comp_index(species_name, pid_name)

Get the component index for a given particle attribute. This is useful to get the corrent ordering of attributes when adding new particles using add_particles().

Parameters
  • species_name (str) – The name of the species

  • pid_name (str) – Name of the component for which the index will be returned

Returns

Integer corresponding to the index of the requested attribute

Return type

int

New components can be added via Python.

pywarpx.picmi.Simulation.extension.add_real_comp(species_name, pid_name, comm=True)

Add a real component to the particle data array.

Parameters
  • species_name (str) – The species name for which the new component will be added

  • pid_name (str) – Name that can be used to identify the new component

  • comm (bool) – Should the component be communicated

Various diagnostics are also accessible from Python. This includes getting the deposited or total charge density from a given species as well as accessing the scraped particle buffer. See the example in Examples/Modules/ParticleBoudaryScrape for a reference on how to interact with scraped particle data.

pywarpx.picmi.Simulation.extension.get_species_charge_sum(species_name, local=False)

Returns the total charge in the simulation due to the given species.

Parameters
  • species_name (str) – The species name for which charge will be summed

  • local (bool) – If True return total charge per processor

pywarpx.picmi.Simulation.extension.depositChargeDensity(species_name, level, clear_rho=True, sync_rho=True)

Deposit the specified species’ charge density in rho_fp in order to access that data via pywarpx.fields.RhoFPWrapper().

Parameters
  • species_name (str) – The species name that will be deposited.

  • level (int) – Which AMR level to retrieve scraped particle data from.

  • clear_rho (bool) – If True, zero out rho_fp before deposition.

  • sync_rho (bool) – If True, perform MPI exchange and properly set boundary cells for rho_fp.

pywarpx.picmi.Simulation.extension.get_particle_boundary_buffer_size(species_name, boundary)

This returns the number of particles that have been scraped so far in the simulation from the specified boundary and of the specified species.

Parameters
  • species_name (str) – Return the number of scraped particles of this species

  • boundary (str) – The boundary from which to get the scraped particle data in the form x/y/z_hi/lo

pywarpx.picmi.Simulation.extension.get_particle_boundary_buffer_size(species_name, boundary)

This returns the number of particles that have been scraped so far in the simulation from the specified boundary and of the specified species.

Parameters
  • species_name (str) – Return the number of scraped particles of this species

  • boundary (str) – The boundary from which to get the scraped particle data in the form x/y/z_hi/lo

pywarpx.picmi.Simulation.extension.get_particle_boundary_buffer_structs(species_name, boundary, level)

This returns a list of numpy arrays containing the particle struct data for a species that has been scraped by a specific simulation boundary. The particle data is represented as a structured numpy array and contains the particle ‘x’, ‘y’, ‘z’, and ‘idcpu’.

The data for the numpy arrays are not copied, but share the underlying memory buffer with WarpX. The numpy arrays are fully writeable.

Parameters
  • species_name (str) – The species name that the data will be returned for

  • boundary (str) – The boundary from which to get the scraped particle data in the form x/y/z_hi/lo or eb.

  • level (int) – Which AMR level to retrieve scraped particle data from.

pywarpx.picmi.Simulation.extension.get_particle_boundary_buffer(species_name, boundary, comp_name, level)

This returns a list of numpy arrays containing the particle array data for a species that has been scraped by a specific simulation boundary.

The data for the numpy arrays are not copied, but share the underlying memory buffer with WarpX. The numpy arrays are fully writeable.

Parameters
  • species_name (str) – The species name that the data will be returned for.

  • boundary (str) – The boundary from which to get the scraped particle data in the form x/y/z_hi/lo or eb.

  • comp_name (str) – The component of the array data that will be returned. If “step_scraped” the special attribute holding the timestep at which a particle was scraped will be returned.

  • level (int) – Which AMR level to retrieve scraped particle data from.

pywarpx.picmi.Simulation.extension.clearParticleBoundaryBuffer()

Clear the buffer that holds the particles lost at the boundaries.

Using Python input as a preprocessor

In this case, only the pure Python version needs to be installed, as described here.

In order to run a new simulation:

  • Create a new directory, where the simulation will be run.

  • Add a Python script in the directory.

The input file should have the line like sim.write_input_file(file_name = 'inputs_from_PICMI') which runs the preprocessor, generating the AMReX inputs file.

  • Run the script with Python:

python <python_script>

where <python_script> is the name of the script. This creates the WarpX input file that you can run as normal with the WarpX executable.