boris_push
- plasmapy.simulation.particle_integrators.boris_push(x, v, B, E, q, m, dt, inplace: bool = True)[source]
The explicit Boris pusher.
Attention
This functionality is under development. Backward incompatible changes might occur in future releases.
- Parameters:
x (
ndarray
) – Particle position at full timestep, in SI (meter) units.v (
ndarray
) – Particle velocity at half timestep, in SI (meter/second) units.B (
ndarray
) – Magnetic field at full timestep, in SI (tesla) units.E (
float
) – Electric field at full timestep, in SI (V/m) units.q (
float
) – Particle charge, in SI (coulomb) units.m (
float
) – Particle mass, in SI (kg) units.dt (
float
) – Timestep, in SI (second) units.inplace (
bool
) – If set, data in the v and x arrays is changed during execution. Defaults to True for performance reasons. If False, this function returnsNone
.
- Returns:
Examples
>>> B = np.array([[0.0, 0.0, 5.0]]) >>> E = np.array([[0.0, 0.0, 0.0]]) >>> x_t0 = np.array([[0.0, 0.0, 0.0]]) >>> v_t0 = np.array([[5.0, 0.0, 0.0]]) >>> x_t1, v_t1 = boris_push( ... x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=False ... ) >>> x_t1 array([[ 0.04993754, -0.00249844, 0. ]]) >>> v_t1 array([[ 4.9937539 , -0.24984385, 0. ]]) >>> boris_push(x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=True) >>> x_t0 array([[ 0.04993754, -0.00249844, 0. ]]) >>> v_t0 array([[ 4.9937539 , -0.24984385, 0. ]]) >>> # B parallel to v >>> B = np.array([[0.0, 0.0, 5.0]]) >>> v_t0 = np.array([[0.0, 0.0, 5.0]]) >>> x_t0 = np.array([[0.0, 0.0, 0.0]]) >>> boris_push(x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=True) >>> # no rotation of vector v >>> v_t0 array([[0., 0., 5.]]) >>> x_t0 array([[0. , 0. , 0.05]]) >>> # B perpendicular to v >>> B = np.array([[5.0, 0.0, 0.0]]) >>> v_t0 = np.array([[0.0, 5.0, 0.0]]) >>> x_t0 = np.array([[0.0, 0.0, 0.0]]) >>> boris_push(x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=True) >>> # rotation of vector v >>> v_t0 array([[ 0. , 4.9937539 , -0.24984385]]) >>> x_t0 array([[ 0. , 0.04993754, -0.00249844]]) >>> # nonzero E and zero B >>> E = np.array([[1.0, 1.0, 1.0]]) >>> B = np.array([[0.0, 0.0, 0.0]]) >>> v_t0 = np.array([[0.0, 0.0, 5.0]]) >>> x_t0 = np.array([[0.0, 0.0, 0.0]]) >>> boris_push(x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=True) >>> v_t0 array([[0.01, 0.01, 5.01]]) >>> x_t0 array([[0.0001, 0.0001, 0.0501]]) >>> boris_push(x=x_t0, v=v_t0, B=B, E=E, q=1.0, m=1.0, dt=0.01, inplace=True) >>> v_t0 array([[0.02, 0.02, 5.02]]) >>> x_t0 array([[0.0003, 0.0003, 0.1003]])
Notes
The Boris algorithm [Boris, 1970] is the standard energy conserving algorithm for particle movement in plasma physics. See Birdsall and Langdon [2004] for more details, and this page on the Boris method for a nice overview.
Conceptually, the algorithm has three phases:
Add half the impulse from electric field.
Rotate the particle velocity about the direction of the magnetic field.
Add the second half of the impulse from the electric field.
This ends up causing the magnetic field action to be properly “centered” in time, and the algorithm, being a symplectic integrator, conserves energy.