XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
push_physics.hpp File Reference
#include "push_controls.hpp"
#include "electric_field.hpp"
#include "perturbed_B_field.hpp"
#include "field_aligned_local_fields.hpp"
#include "push_physics.tpp"
Include dependency graph for push_physics.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

template<class Device , KinType KT, PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void derivs_single_with_e_elec (const Grid< Device > &grid, const PushControls &push_controls, const Species< Device > &species, const MagneticField< Device > &magnetic_field, const GridFieldPack< Device, PIT > &gfpack, const PerturbedBField< Device > &perturbed_B_field, SimdParticles &part, SimdPhase &dy, SimdGridWeights< Order::One, PIT > &grid_wts, const FieldAlignedLocalFields< KT, PIT > &E_mag, double time, Simd< double > *vf_diag=nullptr)
 
template<class Device >
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 SimdVector &efield, SimdPhase &yprime, Simd< double > *vf_diag)
 

Function Documentation

template<class Device >
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 SimdVector efield,
SimdPhase yprime,
Simd< double > *  vf_diag 
)

Get the phase derivatives of a vector of particles

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]efieldVector of electric fields
[in]efield00Vector of 00 electric fields for DELTAF_CONV
[in]ddpotdtVector of dPotential/dt for DELTAF_CONV
[in]f0denVector of f0 density for DELTAF_CONV
[in]f0ddenVector of d[f0 density]/dpsi for DELTAF_CONV
[in]f0tevVector of f0 temperature in eV for DELTAF_CONV
[in]f0dtevVector of d[f0 tempearture]/dpsi for DELTAF_CONV
[in]dAh,Ah,dAs,As: Vector of vector potentials for EXPLICIT_EM. Minus sign included in dAh and dAs
[out]yprimeVector of derivatives of the particle phases
[out]vf_diagOptional vector of quantities for diagnostic

Get the phase derivatives of a vector of particles

Parameters
[in]push_controls
[in]species
[in]partVector of particles
[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]efieldVector of electric fields
[in]dEr_B2Vector of dEr_B2
[in]dEz_B2Vector of dEr_B2
[in]du2_EVector of du2_E
[out]yprimeVector of derivatives of the particle phases
[out]vf_diagOptional vector of quantities for diagnostic

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device , KinType KT, PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void derivs_single_with_e_elec ( const Grid< Device > &  grid,
const PushControls push_controls,
const Species< Device > &  species,
const MagneticField< Device > &  magnetic_field,
const GridFieldPack< Device, PIT > &  gfpack,
const PerturbedBField< Device > &  perturbed_B_field,
SimdParticles part,
SimdPhase dy,
SimdGridWeights< Order::One, PIT > &  grid_wts,
const FieldAlignedLocalFields< KT, PIT > &  E_mag,
double  time,
Simd< double > *  vf_diag 
)

Compute field quantities needed for the push, then output the phase derivatives of a vector of particles

Parameters
[in,out]partVector of particles
[in]dyVector of derivatives of the particle phases
[in]grid_wtsVector of previous triangle index (guess it hasnt moved, to reduce search time)
[in]timeGlobal simulation time, not used at the moment
[out]vf_diagOptional vector of quantities for diagnostic

Here is the call graph for this function:

Here is the caller graph for this function: