ParticleTracker

class plasmapy.simulation.particle_tracker.particle_tracker.ParticleTracker(
grids: AbstractGrid | Iterable[AbstractGrid],
termination_condition: AbstractTerminationCondition | None = None,
save_routine: AbstractSaveRoutine | None = None,
particle_integrator: type[AbstractIntegrator] | None = None,
dt=None,
dt_range=None,
field_weighting: str = 'volume averaged',
req_quantities=None,
verbose: bool = True,
)[source]

Bases: object

A particle tracker for particles in electric and magnetic fields without inter-particle interactions.

Particles are instantiated and pushed through a grid of provided E and B fields using the Boris push algorithm. These fields are specified as part of a grid which are then interpolated to determine the local field acting on each particle.

The time step used in the push routine can be specified, or an adaptive time step will be determined based off the gyroradius of the particle. Some save routines involve time stamping the location and velocities of particles at fixed intervals. In order for this data to be coherent, it is required that all particles follow the same time step. This is referred to as a synchronized time step. If no time step is specified and the provided save routine does not require a synchronized time step, then an adaptive time step is calculated independently for each particle.

The simulation will push particles through the provided fields until a condition is met. This termination condition is provided as an instance of the AbstractTerminationCondition class as arguments to the simulation constructor. The results of a simulation can be exported by specifying an instance of the AbstractSaveRoutine class to the run method.

Parameters:
  • grids (An instance of AbstractGrid) – A Grid object or list of grid objects containing the required quantities. The list of required quantities varies depending on other keywords.

  • termination_condition (AbstractTerminationCondition) – An subclass of AbstractTerminationCondition which determines when the simulation has finished. See AbstractTerminationCondition for more details.

  • save_routine (AbstractSaveRoutine, optional) – An subclass of AbstractSaveRoutine which determines which time steps of the simulation to save. The default is DoNotSaveSaveRoutine. See AbstractSaveRoutine for more details.

  • particle_integrator (AbstractIntegrator, optional) – An subclass of AbstractIntegrator that is responsible for implementing the push behavior of the simulation when provided the electric and magnetic fields. The default value is set to RelativisticBorisIntegrator. See AbstractIntegrator for more information on how to implement custom push routines.

  • dt (Quantity, optional) – An explicitly set time step in units convertible to seconds. Setting this optional keyword overrules the adaptive time step capability and forces the use of this time step throughout.

  • dt_range (tuple of shape (2,) of Quantity, optional) – If specified, the calculated adaptive time step will be clamped between the first and second values.

  • field_weighting (str) –

    String that selects the field weighting algorithm used to determine what fields are felt by the particles. Options are:

    • ’nearest neighbor’: Particles are assigned the fields on

      the grid vertex closest to them.

    • ’volume averaged’The fields experienced by a particle are a

      volume-average of the eight grid points surrounding them.

    The default is ‘volume averaged’.

  • req_quantities (list of str, default : None) – A list of quantity keys required to be specified on the Grid object. The base particle pushing simulation requires the quantities [E_x, E_y, E_z, B_x, B_y, B_z]. This keyword is for specifying quantities in addition to these six. If any required quantities are missing, those quantities will be assumed to be zero everywhere. A warning will be raised if any of the additional required quantities are missing and are set to zero.

  • verbose (bool, optional) – If true, updates on the status of the program will be printed into the standard output while running. The default is True.

Warns:

`~plasmapy.utils.exceptions.RelativityWarning` – The provided integrator does not account for relativistic corrections and particles have reached a non-negligible fraction of the speed of light.

Notes

We adopt the convention of NaN values to represent various states of a particle.

If the particle’s position and velocity are not NaN, the particle is being tracked and evolved. If the particle’s position is not NaN, but the velocity is NaN, the particle has been stopped (i.e. it is still in the simulation but is no longer evolved.) If both the particle’s position and velocity are set to NaN, then the particle has been removed from the simulation.

Example

>>> from plasmapy.particles import Particle
>>> from plasmapy.plasma.grids import CartesianGrid
>>> from plasmapy.simulation.particle_tracker.particle_tracker import (
...     ParticleTracker,
... )
>>> from plasmapy.simulation.particle_tracker.save_routines import (
...     IntervalSaveRoutine,
... )
>>> from plasmapy.simulation.particle_tracker.termination_conditions import (
...     TimeElapsedTerminationCondition,
... )
>>> import astropy.units as u
>>> import numpy as np
>>> example_particle = Particle("p+")
>>> grid = CartesianGrid(-1e6 * u.m, 1e6 * u.m, num=2)
>>> grid_shape = (2, 2, 2)
>>> Bz = np.full(grid_shape, 1) * u.T
>>> grid.add_quantities(B_z=Bz)
>>> x = [[0, 0, 0]] * u.m
>>> v = [[1, 0, 0]] * u.m / u.s
>>> termination_condition = TimeElapsedTerminationCondition(6.28 * u.s)
>>> save_routine = IntervalSaveRoutine(6.28 * u.s / 10)
>>> simulation = ParticleTracker(
...     grid,
...     termination_condition,
...     dt=1e-2 * u.s,
...     save_routine=save_routine,
...     field_weighting="nearest neighbor",
... )
>>> simulation.load_particles(x, v, example_particle)
>>> simulation.run()
>>> print(
...     simulation.save_routine.results["time"][-1],
...     simulation.save_routine.results["x"][-1],
... )
6.29999999999991 s [[-1.73302335e-08  1.31539878e-05  0.00000000e+00]] m

Attributes Summary

fract_entered

The fraction of particles that have entered the grid.

is_adaptive_time_step

Return whether the simulation is calculating an adaptive time step or using the user-provided time step.

is_synchronized_time_step

Return if the simulation is applying the same time step across all particles.

nparticles_removed

Return the number of particles currently being tracked.

nparticles_stopped

Return the number of particles currently being tracked.

nparticles_tracked

Return the number of particles currently being tracked.

num_entered

Count the number of particles that have entered the grids.

num_grids

The number of grids specified at instantiation.

on_any_grid

Binary array for each particle indicating whether it is currently on ANY grid.

particles_on_grid

Returns a boolean mask of shape [ngrids, nparticles] corresponding to whether or not the particle is on the associated grid.

vmax

The maximum velocity of any particle in the simulation.

Methods Summary

add_stopping(method[, materials, I])

Enable particle stopping using experimental stopping powers.

load_particles(x, v, particle)

Load arrays of particle positions and velocities.

run()

Runs a particle-tracing simulation.

setup_adaptive_time_step([...])

Set parameters for the adaptive time step candidates.

Attributes Documentation

fract_entered

The fraction of particles that have entered the grid. The denominator of this fraction is based off the number of tracked particles, and therefore does not include stopped or removed particles.

is_adaptive_time_step

Return whether the simulation is calculating an adaptive time step or using the user-provided time step.

is_synchronized_time_step

Return if the simulation is applying the same time step across all particles.

nparticles_removed

Return the number of particles currently being tracked. That is, they do not have NaN position or velocity.

nparticles_stopped

Return the number of particles currently being tracked. That is, they do not have NaN position or velocity.

nparticles_tracked

Return the number of particles currently being tracked. That is, they do not have NaN position or velocity.

num_entered

Count the number of particles that have entered the grids. This number is calculated by summing the number of non-zero entries in the entered grid array.

num_grids

The number of grids specified at instantiation.

on_any_grid

Binary array for each particle indicating whether it is currently on ANY grid.

particles_on_grid

Returns a boolean mask of shape [ngrids, nparticles] corresponding to whether or not the particle is on the associated grid.

vmax

The maximum velocity of any particle in the simulation.

This quantity is used for determining the grid crossing maximum time step.

Methods Documentation

add_stopping(
method: Literal['NIST', 'Bethe'],
materials: list[str] | None = None,
I: Quantity | None = None,
)[source]

Enable particle stopping using experimental stopping powers.

Interpolation of stopping powers is conducted using data from the NIST PSTAR database. This information is combined with the mass density quantity provided in the grids to calculate the energy loss over the distance travelled during a timestep.

load_particles(
x,
v,
particle: Particle,
) None[source]

Load arrays of particle positions and velocities.

Parameters:
  • x (Quantity, shape (N,3)) – Positions for N particles

  • v (Quantity, shape (N,3)) – Velocities for N particles

  • particle (particle-like) – Representation of the particle species as either a Particle object or a string representation.

run() None[source]

Runs a particle-tracing simulation. Time steps are adaptively calculated based on the local grid resolution of the particles and the electric and magnetic fields they are experiencing.

Return type:

None

setup_adaptive_time_step(
time_steps_per_gyroperiod: int | None = 12,
Courant_parameter: float | None = 0.5,
) None[source]

Set parameters for the adaptive time step candidates.

Parameters:
  • time_steps_per_gyroperiod (int, optional) – The minimum ratio of the particle gyroperiod to the timestep. Higher numbers correspond to higher temporal resolution. The default is twelve.

  • Courant_parameter (float, optional) – The Courant parameter is the minimum ratio of the timestep to the grid crossing time, grid cell length / particle velocity. Lower Courant numbers correspond to higher temporal resolution.

Notes

Two candidates are calculated for the adaptive time step: a time step based on the gyroradius of the particle and a time step based on the resolution of the grid. The candidate associated with the gyroradius of the particle takes a time_steps_per_gyroperiod parameter that specifies how many times the orbit of a gyrating particles will be subdivided. The other candidate, associated with the spatial resolution of the grid object, calculates a time step using the time it would take the fastest particle to cross some fraction of a grid cell length. This fraction is the Courant number.