XGC1
|
#include <perturbed_B_field.hpp>
Public Member Functions | |
PerturbedBField (NLReader::NamelistReader &nlr, bool rampup_vac_in, int num_ntor_in, int rampup_time_in, int start_time_in, int nnode_in, int *ntor_in) | |
PerturbedBField () | |
void | update_rampup_fac_from_fortran () const |
void | copy_bfield_from_fortran () const |
KOKKOS_INLINE_FUNCTION void | get_delta_b (const Grid< Device > &grid, const Simd< double > &phi, const SimdGridWeights< Order::One, PhiInterpType::None > &grid_wts, SimdVector &tdb) const |
Public Attributes | |
bool | rampup_vac |
(.true.) Ramp up perturbed field slowly, (.false.) turn on perturbed field abruptly More... | |
int | num_ntor |
Number of toroidal mode numbers (<= sml_nphi_total) More... | |
int | rampup_time |
Number of time steps over which the perturbed field is ramped up. More... | |
int | start_time |
Time step in which perturbed field is switched on. More... | |
int | mode |
Mode of RMP operation (details in XGC_core/module_ptb_3db.F90 –> ptb_3db_mode) More... | |
bool | full_spec_on |
(EM-only) whether to retain the full RMP spectrum from M3D-C1 or only |m/q-n|<=sml_mode_select_mres_q More... | |
double | rampup_fac |
The current relative amplitude of the vacuum vector potential (if ptb_3db_mode==2) - this is updated during the simulation. More... | |
Kokkos::View< int *, Kokkos::LayoutRight, Device > | ntor |
Array to store the toroidal mode numbers. More... | |
Kokkos::View< Field < VarType::Vector, PhiInterpType::None > **, Kokkos::LayoutRight, Device > | bfield_im_vac |
perturbed vacuum field on XGC mesh, imaginary part More... | |
Kokkos::View< Field < VarType::Vector, PhiInterpType::None > **, Kokkos::LayoutRight, Device > | bfield_re_vac |
perturbed vacuum field on XGC mesh, real part, dimensions: (grid vertex, R-Z-phi components, tor. mode number) More... | |
PerturbedBField< Device >::PerturbedBField | ( | NLReader::NamelistReader & | nlr, |
bool | rampup_vac_in, | ||
int | num_ntor_in, | ||
int | rampup_time_in, | ||
int | start_time_in, | ||
int | nnode_in, | ||
int * | ntor_in | ||
) |
Constructor for perturbed magnetic field class
|
inline |
|
inline |
KOKKOS_INLINE_FUNCTION void PerturbedBField< Device >::get_delta_b | ( | const Grid< Device > & | grid, |
const Simd< double > & | phi, | ||
const SimdGridWeights< Order::One, PhiInterpType::None > & | grid_wts, | ||
SimdVector & | tdb | ||
) | const |
Get magnetic field perturbation vector for a vector of particle locations
[in] | grid | The grid object |
[in] | phi | Phi coordinate |
[in] | grid_wts | Triangle the particles are in |
[out] | tdb | Magnetic field perturbation vector |
|
inline |
Kokkos::View<Field<VarType::Vector,PhiInterpType::None>**,Kokkos::LayoutRight,Device> PerturbedBField< Device >::bfield_im_vac |
perturbed vacuum field on XGC mesh, imaginary part
Kokkos::View<Field<VarType::Vector,PhiInterpType::None>**,Kokkos::LayoutRight,Device> PerturbedBField< Device >::bfield_re_vac |
perturbed vacuum field on XGC mesh, real part, dimensions: (grid vertex, R-Z-phi components, tor. mode number)
bool PerturbedBField< Device >::full_spec_on |
(EM-only) whether to retain the full RMP spectrum from M3D-C1 or only |m/q-n|<=sml_mode_select_mres_q
int PerturbedBField< Device >::mode |
Mode of RMP operation (details in XGC_core/module_ptb_3db.F90 –> ptb_3db_mode)
Kokkos::View<int*,Kokkos::LayoutRight,Device> PerturbedBField< Device >::ntor |
Array to store the toroidal mode numbers.
int PerturbedBField< Device >::num_ntor |
Number of toroidal mode numbers (<= sml_nphi_total)
double PerturbedBField< Device >::rampup_fac |
The current relative amplitude of the vacuum vector potential (if ptb_3db_mode==2) - this is updated during the simulation.
int PerturbedBField< Device >::rampup_time |
Number of time steps over which the perturbed field is ramped up.
bool PerturbedBField< Device >::rampup_vac |
(.true.) Ramp up perturbed field slowly, (.false.) turn on perturbed field abruptly
int PerturbedBField< Device >::start_time |
Time step in which perturbed field is switched on.