WarpX
Functions
PushSelector.H File Reference
#include "Particles/Pusher/UpdateMomentumBoris.H"
#include "Particles/Pusher/UpdateMomentumBorisWithRadiationReaction.H"
#include "Particles/Pusher/UpdateMomentumHigueraCary.H"
#include "Particles/Pusher/UpdateMomentumVay.H"
#include "Particles/WarpXParticleContainer.H"
#include "Utils/WarpXAlgorithmSelection.H"
#include <AMReX_REAL.H>
#include <limits>

Go to the source code of this file.

Functions

template<int do_sync>
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void doParticleMomentumPush (amrex::ParticleReal &ux, amrex::ParticleReal &uy, amrex::ParticleReal &uz, const amrex::ParticleReal Ex, const amrex::ParticleReal Ey, const amrex::ParticleReal Ez, const amrex::ParticleReal Bx, const amrex::ParticleReal By, const amrex::ParticleReal Bz, const int ion_lev, const amrex::ParticleReal m, const amrex::ParticleReal a_q, const int pusher_algo, const int do_crr, const amrex::Real t_chi_max, const amrex::Real dt)
 Push momentum for a single particle. More...
 

Function Documentation

◆ doParticleMomentumPush()

template<int do_sync>
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void doParticleMomentumPush ( amrex::ParticleReal &  ux,
amrex::ParticleReal &  uy,
amrex::ParticleReal &  uz,
const amrex::ParticleReal  Ex,
const amrex::ParticleReal  Ey,
const amrex::ParticleReal  Ez,
const amrex::ParticleReal  Bx,
const amrex::ParticleReal  By,
const amrex::ParticleReal  Bz,
const int  ion_lev,
const amrex::ParticleReal  m,
const amrex::ParticleReal  a_q,
const int  pusher_algo,
const int  do_crr,
const amrex::Real  t_chi_max,
const amrex::Real  dt 
)

Push momentum for a single particle.

Template Parameters
do_syncWhether to include quantum synchrotron radiation (QSR)
Parameters
ux,uy,uzParticle momentum
Ex,Ey,EzElectric field on particles.
Bx,By,BzMagnetic field on particles.
ion_levIonization level of this particle (0 if ionization not on)
mMass of this species.
a_qCharge of this species.
pusher_algo0: Boris, 1: Vay, 2: HigueraCary
do_crrWhether to do the classical radiation reaction
t_chi_maxCutoff chi for QSR
dtTime step size