8 #ifndef WARPX_PARTICLES_PUSHER_UPDATEPOSITION_H_
9 #define WARPX_PARTICLES_PUSHER_UPDATEPOSITION_H_
24 void UpdatePosition(amrex::ParticleReal& x, amrex::ParticleReal& y, amrex::ParticleReal& z,
25 const amrex::ParticleReal ux,
const amrex::ParticleReal uy,
const amrex::ParticleReal uz,
26 const amrex::Real
dt )
28 using namespace amrex::literals;
33 const amrex::ParticleReal inv_gamma = 1._prt/std::sqrt(1._prt + (ux*ux + uy*uy + uz*uz)*inv_c2);
35 #if (AMREX_SPACEDIM >= 2)
36 x += ux * inv_gamma *
dt;
40 #if defined(WARPX_DIM_3D) || defined(WARPX_DIM_RZ)
41 y += uy * inv_gamma *
dt;
45 z += uz * inv_gamma *
dt;
56 const amrex::ParticleReal ux_n,
const amrex::ParticleReal uy_n,
const amrex::ParticleReal uz_n,
57 const amrex::ParticleReal ux,
const amrex::ParticleReal uy,
const amrex::ParticleReal uz,
58 const amrex::Real
dt )
60 using namespace amrex::literals;
66 const amrex::ParticleReal ux_np1 = 2._prt*ux - ux_n;
67 const amrex::ParticleReal uy_np1 = 2._prt*uy - uy_n;
68 const amrex::ParticleReal uz_np1 = 2._prt*uz - uz_n;
69 const amrex::ParticleReal gamma_n = std::sqrt(1._prt + (ux_n*ux_n + uy_n*uy_n + uz_n*uz_n)*inv_c2);
70 const amrex::ParticleReal gamma_np1 = std::sqrt(1._prt + (ux_np1*ux_np1 + uy_np1*uy_np1 + uz_np1*uz_np1)*inv_c2);
71 const amrex::ParticleReal inv_gamma = 2.0_prt/(gamma_n + gamma_np1);
74 #if (AMREX_SPACEDIM >= 2)
75 x += ux * inv_gamma *
dt;
79 #if defined(WARPX_DIM_3D) || defined(WARPX_DIM_RZ)
80 y += uy * inv_gamma *
dt;
84 z += uz * inv_gamma *
dt;
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_INLINE void UpdatePositionImplicit(amrex::ParticleReal &x, amrex::ParticleReal &y, amrex::ParticleReal &z, const amrex::ParticleReal ux_n, const amrex::ParticleReal uy_n, const amrex::ParticleReal uz_n, const amrex::ParticleReal ux, const amrex::ParticleReal uy, const amrex::ParticleReal uz, const amrex::Real dt)
Push the particle's positions over one timestep, given the value of its momenta ux,...
Definition: UpdatePosition.H:55
AMREX_GPU_HOST_DEVICE AMREX_INLINE void UpdatePosition(amrex::ParticleReal &x, amrex::ParticleReal &y, amrex::ParticleReal &z, const amrex::ParticleReal ux, const amrex::ParticleReal uy, const amrex::ParticleReal uz, const amrex::Real dt)
Push the particle's positions over one timestep, given the value of its momenta ux,...
Definition: UpdatePosition.H:24
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
float dt
Definition: stencil.py:442