XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | List of all members
PerturbedBField< Device > Class Template Reference

#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
 
KOKKOS_INLINE_FUNCTION void get_delta_b (const Grid< Device > &grid, const Simd< double > &fld_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...
 

Constructor & Destructor Documentation

template<class Device >
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

Here is the call graph for this function:

template<class Device>
PerturbedBField< Device >::PerturbedBField ( )
inline

Member Function Documentation

template<class Device >
KOKKOS_INLINE_FUNCTION void PerturbedBField< Device >::get_delta_b ( const Grid< Device > &  grid,
const Simd< double > &  fld_phi,
const SimdGridWeights< Order::One, PhiInterpType::None > &  grid_wts,
SimdVector tdb 
) const

Get magnetic field perturbation vector for a vector of particle locations

Parameters
[in]gridThe grid object
[in]fld_phiPhi coordinate
[in]grid_wtsTriangle the particles are in
[out]tdbMagnetic field perturbation vector

Here is the caller graph for this function:

template<class Device>
void PerturbedBField< Device >::update_rampup_fac_from_fortran ( ) const
inline

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

template<class Device>
Kokkos::View<Field<VarType::Vector,PhiInterpType::None>**,Kokkos::LayoutRight,Device> PerturbedBField< Device >::bfield_im_vac

perturbed vacuum field on XGC mesh, imaginary part

template<class Device>
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)

template<class Device>
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

template<class Device>
int PerturbedBField< Device >::mode

Mode of RMP operation (details in XGC_core/module_ptb_3db.F90 –> ptb_3db_mode)

template<class Device>
Kokkos::View<int*,Kokkos::LayoutRight,Device> PerturbedBField< Device >::ntor

Array to store the toroidal mode numbers.

template<class Device>
int PerturbedBField< Device >::num_ntor

Number of toroidal mode numbers (<= sml_nphi_total)

template<class Device>
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.

template<class Device>
int PerturbedBField< Device >::rampup_time

Number of time steps over which the perturbed field is ramped up.

template<class Device>
bool PerturbedBField< Device >::rampup_vac

(.true.) Ramp up perturbed field slowly, (.false.) turn on perturbed field abruptly

template<class Device>
int PerturbedBField< Device >::start_time

Time step in which perturbed field is switched on.


The documentation for this class was generated from the following files: