XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
em_push_physics.tpp File Reference
#include "globals.hpp"
#include "push_controls.hpp"
#include "magnetic_field.hpp"
#include "grid.hpp"
#include "particles.hpp"
#include "species.hpp"
#include "perturbed_B_field.hpp"
Include dependency graph for em_push_physics.tpp:

Functions

KOKKOS_INLINE_FUNCTION void get_ExB_drift (double D, double drift_on, const SimdVector &bfield, const SimdVector &E, int i_simd, double over_B2, double inv_r, double(&yp_exb)[3])
 
KOKKOS_INLINE_FUNCTION double get_denergy_dt_es (double charge, int i_simd, const SimdVector &bfield, const SimdVector &tdb, const SimdVector &E, double rho, double vp, double mu, double B, double c_m, double over_B, double over_B2, double dbdr, double dbdz, double dbdphi, const double(&curl_B)[3], const double drift_on, double D)
 
KOKKOS_INLINE_FUNCTION double get_denergy_dt_em (double mass, int i_simd, const Simd< double > &f0flow, const Simd< double > &Epar_em, const double(&Bs)[3], const double(&Bsa)[3], const double(&dXdt)[3], double fr_pert, double fz_pert, double fp_pert, double fr_mag, double fz_mag, double fp_mag, double vp, double mu, double B, double c_m, double over_B, double dbdr, double dbdz, double dbdphi, double D)
 
template<MarkerType MT>
KOKKOS_INLINE_FUNCTION double get_dfdp (double mass, double vp, double mu, double B, double inv_temp, const LocalEquilProfiles< MT > &local_equil_profiles, int i_simd)
 
KOKKOS_INLINE_FUNCTION double get_dweight (double dfdp, const SimdVector2D &gradpsi, int i_simd, const double(&yp)[3], double dEdtT)
 
template<MarkerType MT>
KOKKOS_INLINE_FUNCTION double get_dw1 (double mass, double charge, bool is_electron, const PushControls &push_controls, int i_simd, double w2, double r, const SimdVector &bfield, const SimdVector2D &gradpsi, const SimdVector &tdb, const LocalFields &fld, const LocalEquilProfiles< MT > &local_equil_profiles, const double(&Bs)[3], const double(&Bsa)[3], double fr_pert, double fz_pert, double fp_pert, double fr2, double fz2, double fp2, double inv_r, double rho, double vp, double mu, double B, double c_m, double over_B, double over_B2, double dbdr, double dbdz, double dbdphi, const double(&curl_B)[3], const double drift_on, double D, double yp_r, double yp_z, double yp_p)
 
KOKKOS_INLINE_FUNCTION void get_yp_pert (double D, const double(&Bs)[3], const double(&Bsa)[3], const SimdVector &bfield, double vp, double c_m, const Simd< double > &Ah, double fr_pert, double fz_pert, double fp_pert, double over_B, double over_B2, double inv_r, int i_simd, double(&yp_pert)[3])
 
KOKKOS_INLINE_FUNCTION double get_adiabatic_response (const SimdVector &E, const SimdVector &E00, const Simd< double > &ddpotdt, double yp_r, double yp_z, double yp_p, double r, double inv_temp, int i_simd)
 
template<>
KOKKOS_INLINE_FUNCTION double get_dw1< MarkerType::ReducedDeltaF > (double mass, double charge, bool is_electron, const PushControls &push_controls, int i_simd, double w2, double r, const SimdVector &bfield, const SimdVector2D &gradpsi, const SimdVector &tdb, const LocalFields &fld, const LocalEquilProfiles< MarkerType::ReducedDeltaF > &local_equil_profiles, const double(&Bs)[3], const double(&Bsa)[3], double fr_pert, double fz_pert, double fp_pert, double fr2, double fz2, double fp2, double inv_r, double rho, double vp, double mu, double B, double c_m, double over_B, double over_B2, double dbdr, double dbdz, double dbdphi, const double(&curl_B)[3], const double drift_on, double D, double yp_r, double yp_z, double yp_p)
 
KOKKOS_INLINE_FUNCTION double get_nb_curl_nb (double over_B2, const SimdVector &bfield, const SimdVector(&jacb)[3], double inv_r, int i_simd)
 
template<class Device , PushDiagToggle PDT>
KOKKOS_INLINE_FUNCTION void derivs_elec (const PushControls &push_controls, const Species< Device > &species, const SimdParticles &part, const Simd< double > &inv_r_vec, const SimdVector &bfield, const Simd< double > &B_mag, const SimdVector(&jacb)[3], const Simd< double > &psi, const SimdVector2D &gradpsi, const SimdVector &tdb, const LocalFields &fld, const LocalEquilProfiles< MT_GLOBAL > &local_equil_profiles, SimdPhase &yprime, VFDiag< PDT > &vf_diag)
 

Function Documentation

template<class Device , PushDiagToggle PDT>
KOKKOS_INLINE_FUNCTION void derivs_elec ( const PushControls push_controls,
const Species< Device > &  species,
const SimdParticles part,
const Simd< double > &  inv_r_vec,
const SimdVector bfield,
const Simd< double > &  B_mag,
const SimdVector(&)  jacb[3],
const Simd< double > &  psi,
const SimdVector2D gradpsi,
const SimdVector tdb,
const LocalFields fld,
const LocalEquilProfiles< MT_GLOBAL > &  local_equil_profiles,
SimdPhase yprime,
VFDiag< PDT > &  vf_diag 
)

Get the phase derivatives of a vector of particles

Parameters
[in]push_controlsSet of options/parameters used in the particle push
[in]speciesContains species parameters
[in]partVector of particles
[in]inv_r_vecVector of 1/r (1/r_axis in cylindrical limit)
[in]bfieldVector of magnetic fields
[in]B_magVector of bfield magnitudes
[in]jacbVector of jacobians
[in]psiVector of psi coordinates
[in]gradpsiVector of psi gradients (r,z)
[in]tdbVector of perturbed magnetic field
[in]fldContains electromagnetic fields at a set of points
[in]local_equil_profilesContains f0 profile information is required for conventional delta-f method (in dw/dt)
[out]yprimeVector of derivatives of the particle phases
[out]vf_diagOptional vector of quantities for diagnostic

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION double get_adiabatic_response ( const SimdVector E,
const SimdVector E00,
const Simd< double > &  ddpotdt,
double  yp_r,
double  yp_z,
double  yp_p,
double  r,
double  inv_temp,
int  i_simd 
)

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION double get_denergy_dt_em ( double  mass,
int  i_simd,
const Simd< double > &  f0flow,
const Simd< double > &  Epar_em,
const double(&)  Bs[3],
const double(&)  Bsa[3],
const double(&)  dXdt[3],
double  fr_pert,
double  fz_pert,
double  fp_pert,
double  fr_mag,
double  fz_mag,
double  fp_mag,
double  vp,
double  mu,
double  B,
double  c_m,
double  over_B,
double  dbdr,
double  dbdz,
double  dbdphi,
double  D 
)

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION double get_denergy_dt_es ( double  charge,
int  i_simd,
const SimdVector bfield,
const SimdVector tdb,
const SimdVector E,
double  rho,
double  vp,
double  mu,
double  B,
double  c_m,
double  over_B,
double  over_B2,
double  dbdr,
double  dbdz,
double  dbdphi,
const double(&)  curl_B[3],
const double  drift_on,
double  D 
)

Here is the caller graph for this function:

template<MarkerType MT>
KOKKOS_INLINE_FUNCTION double get_dfdp ( double  mass,
double  vp,
double  mu,
double  B,
double  inv_temp,
const LocalEquilProfiles< MT > &  local_equil_profiles,
int  i_simd 
)

Here is the caller graph for this function:

template<MarkerType MT>
KOKKOS_INLINE_FUNCTION double get_dw1 ( double  mass,
double  charge,
bool  is_electron,
const PushControls push_controls,
int  i_simd,
double  w2,
double  r,
const SimdVector bfield,
const SimdVector2D gradpsi,
const SimdVector tdb,
const LocalFields fld,
const LocalEquilProfiles< MT > &  local_equil_profiles,
const double(&)  Bs[3],
const double(&)  Bsa[3],
double  fr_pert,
double  fz_pert,
double  fp_pert,
double  fr2,
double  fz2,
double  fp2,
double  inv_r,
double  rho,
double  vp,
double  mu,
double  B,
double  c_m,
double  over_B,
double  over_B2,
double  dbdr,
double  dbdz,
double  dbdphi,
const double(&)  curl_B[3],
const double  drift_on,
double  D,
double  yp_r,
double  yp_z,
double  yp_p 
)

Here is the caller graph for this function:

template<>
KOKKOS_INLINE_FUNCTION double get_dw1< MarkerType::ReducedDeltaF > ( double  mass,
double  charge,
bool  is_electron,
const PushControls push_controls,
int  i_simd,
double  w2,
double  r,
const SimdVector bfield,
const SimdVector2D gradpsi,
const SimdVector tdb,
const LocalFields fld,
const LocalEquilProfiles< MarkerType::ReducedDeltaF > &  local_equil_profiles,
const double(&)  Bs[3],
const double(&)  Bsa[3],
double  fr_pert,
double  fz_pert,
double  fp_pert,
double  fr2,
double  fz2,
double  fp2,
double  inv_r,
double  rho,
double  vp,
double  mu,
double  B,
double  c_m,
double  over_B,
double  over_B2,
double  dbdr,
double  dbdz,
double  dbdphi,
const double(&)  curl_B[3],
const double  drift_on,
double  D,
double  yp_r,
double  yp_z,
double  yp_p 
)

Here is the call graph for this function:

KOKKOS_INLINE_FUNCTION double get_dweight ( double  dfdp,
const SimdVector2D gradpsi,
int  i_simd,
const double(&)  yp[3],
double  dEdtT 
)

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION void get_ExB_drift ( double  D,
double  drift_on,
const SimdVector bfield,
const SimdVector E,
int  i_simd,
double  over_B2,
double  inv_r,
double(&)  yp_exb[3] 
)

Get the phase derivatives of a vector of particles

Parameters
[in]D
[in]drift_on
[in]bfieldVector of magnetic fields
[in]E
[in]i_simdvector index
[in]over_B21/|B|^2
[in]inv_r1/r
[out]yp_exb

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION double get_nb_curl_nb ( double  over_B2,
const SimdVector bfield,
const SimdVector(&)  jacb[3],
double  inv_r,
int  i_simd 
)

Get normalized b dot curl of normalized b

Parameters
[in]over_B21/|B|^2
[in]bfieldVector of magnetic fields
[in]jacbVector of jacobians
[in]inv_r1/r
[in]i_simdvector index
Returns
double normalized b dot curl of normalized b

Here is the caller graph for this function:

KOKKOS_INLINE_FUNCTION void get_yp_pert ( double  D,
const double(&)  Bs[3],
const double(&)  Bsa[3],
const SimdVector bfield,
double  vp,
double  c_m,
const Simd< double > &  Ah,
double  fr_pert,
double  fz_pert,
double  fp_pert,
double  over_B,
double  over_B2,
double  inv_r,
int  i_simd,
double(&)  yp_pert[3] 
)

Here is the caller graph for this function: