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

#include <magnetic_field.hpp>

Collaboration diagram for MagneticField< Device >:
Collaboration graph
[legend]

Public Member Functions

 MagneticField (double bt_sign_in, double bp_sign_in, int ff_step_in, int ff_order_in, double bd_min_r_in, double bd_max_r_in, double bd_min_z_in, double bd_max_z_in, double inpsi_in, double outpsi_in, double *rc_in, double *zc_in, BicubCoeff *acoeff_in, int nr_in, int nz_in, double rmin_in, double zmin_in, double dr_inv_in, double dz_inv_in, OneDCoeff *one_d_cub_acoef_in, int ncoeff_in, double max_psi_in, double min_psi_in, double one_d_cub_dpsi_inv_in, double eq_min_r, double eq_max_r, double eq_min_z, double eq_max_z, double eq_x_psi, double epsil_psi, double eq_x_r, double eq_x_slope, double eq_x_z, double eq_x2_r, double eq_x2_slope, double eq_x2_z, double eq_x2_psi, double eq_axis_r, double eq_axis_z, int eq_mpsi)
 
 MagneticField (NLReader::NamelistReader &nlr, const MagneticEquilFiles &equil_files, const View< Equil::XPoint *, HostType > &xpts)
 
template<class Device2 >
MagneticField< Device2 > mirror () const
 
 MagneticField (PsiOption psi_opt, double safety_factor_coeff=1.0, int eq_mr_in=-1, int eq_mz_in=-1, int eq_mpsi_in=-1)
 
 MagneticField ()
 
KOKKOS_INLINE_FUNCTION void get_psi (const SimdVector2D &x, Simd< double > &psi_out) const
 
KOKKOS_INLINE_FUNCTION double get_psi (double r, double z) const
 
KOKKOS_INLINE_FUNCTION void get_psi_unit_vec (const SimdVector2D &x, Simd< double > &cosa, Simd< double > &sina) const
 
KOKKOS_INLINE_FUNCTION double I_value (double psi_in, int rgn3) const
 
KOKKOS_INLINE_FUNCTION double I_deriv (double psi_in, int rgn3) const
 
KOKKOS_INLINE_FUNCTION double geometry_r (double r) const
 
KOKKOS_INLINE_FUNCTION void field (const SimdVector2D &x, SimdVector &bvec, SimdVector(&jacb)[3], Simd< double > &psivec, SimdVector2D &gradpsi, SimdVector &tdb, Simd< bool > &rz_outside) const
 
KOKKOS_INLINE_FUNCTION void bvec_interpol (double r, double z, double phi, double &br, double &bz, double &bphi) const
 
KOKKOS_INLINE_FUNCTION void bmag_interpol (const SimdVector &v, Simd< double > &bmag) const
 
KOKKOS_INLINE_FUNCTION void bmag_interpol (const SimdVector2D &x, const Simd< double > &phi, Simd< double > &bmag) const
 
KOKKOS_INLINE_FUNCTION void follow_field (const SimdVector2D &x_org, const Simd< double > &phi_org, const Simd< double > &phi_dest, SimdVector2D &x_dest) const
 
template<>
 MagneticField (NLReader::NamelistReader &nlr, const MagneticEquilFiles &equil_files, const View< Equil::XPoint *, HostType > &xpts)
 

Public Attributes

double bt_sign
 Whether toroidal field is reversed? More...
 
double bp_sign
 Whether poloidal field is reversed? More...
 
int ff_step
 Number of steps taken when projecting the particle location onto the midplane. More...
 
int ff_order
 Order of RK scheme used for field following. Can be 1, 2, or 4. More...
 
RZBounds bounds
 Simulation boundary. More...
 
double inpsi
 Boundary condition used in a few spots. More...
 
double outpsi
 Boundary condition used in a few spots. More...
 
Bicub< Device > psi_bicub
 The object for interpolating psi (magnetic flux surfaces) More...
 
CubInterp< Device > I_interp
 The object for interpolating I (for deriving toroidal magnetic field) More...
 
Equilibrium equil
 The object containing information about the magnetic equilibrium. More...
 

Private Member Functions

KOKKOS_INLINE_FUNCTION void derivs (const double(&x)[2], double phi, double(&dx)[2]) const
 

Constructor & Destructor Documentation

template<class Device >
MagneticField< Device >::MagneticField ( double  bt_sign_in,
double  bp_sign_in,
int  ff_step_in,
int  ff_order_in,
double  bd_min_r_in,
double  bd_max_r_in,
double  bd_min_z_in,
double  bd_max_z_in,
double  inpsi_in,
double  outpsi_in,
double *  rc_in,
double *  zc_in,
BicubCoeff acoeff_in,
int  nr_in,
int  nz_in,
double  rmin_in,
double  zmin_in,
double  dr_inv_in,
double  dz_inv_in,
OneDCoeff one_d_cub_acoef_in,
int  ncoeff_in,
double  max_psi_in,
double  min_psi_in,
double  one_d_cub_dpsi_inv_in,
double  eq_min_r,
double  eq_max_r,
double  eq_min_z,
double  eq_max_z,
double  eq_x_psi,
double  epsil_psi,
double  eq_x_r,
double  eq_x_slope,
double  eq_x_z,
double  eq_x2_r,
double  eq_x2_slope,
double  eq_x2_z,
double  eq_x2_psi,
double  eq_axis_r,
double  eq_axis_z,
int  eq_mpsi 
)

Constructor for magnetic field

template<class Device>
MagneticField< Device >::MagneticField ( NLReader::NamelistReader nlr,
const MagneticEquilFiles equil_files,
const View< Equil::XPoint *, HostType > &  xpts 
)
inline
template<class Device>
MagneticField< Device >::MagneticField ( PsiOption  psi_opt,
double  safety_factor_coeff = 1.0,
int  eq_mr_in = -1,
int  eq_mz_in = -1,
int  eq_mpsi_in = -1 
)
inline

Here is the call graph for this function:

template<class Device>
MagneticField< Device >::MagneticField ( )
inline
template<>
MagneticField< HostType >::MagneticField ( NLReader::NamelistReader nlr,
const MagneticEquilFiles equil_files,
const View< Equil::XPoint *, HostType > &  xpts 
)
inline

Here is the call graph for this function:

Member Function Documentation

template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::bmag_interpol ( const SimdVector v,
Simd< double > &  bmag 
) const

For vector of coordinates, get the B-field magnitude

Parameters
[in]vVector of {r,z,phi} coordinates
[out]bmagVector of magnetic field magnitudes

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::bmag_interpol ( const SimdVector2D x,
const Simd< double > &  phi,
Simd< double > &  bmag 
) const

For vector of coordinates, get the B-field magnitude

Parameters
[in]xVector of {r,z} coordinates
[out]bmagVector of magnetic field magnitudes
template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::bvec_interpol ( double  r,
double  z,
double  phi,
double &  br,
double &  bz,
double &  bphi 
) const

Get the magnetic field components at single location

Parameters
[in]rr coordinate
[in]zz coordinate
[in]phiphi coordinate
[out]brr component of B
[out]bzz component of B
[out]bphiphi component of B

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::derivs ( const double(&)  x[2],
double  phi,
double(&)  dx[2] 
) const
private

Get the derivatives of the magnetic field to see where the particle should move to as it traces the field

Parameters
[in]x(r,z) coordinates
[in]phiphi coordinate
[out]dxchange in (r,z)
template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::field ( const SimdVector2D x,
SimdVector bvec,
SimdVector(&)  jacb[3],
Simd< double > &  psivec,
SimdVector2D gradpsi,
SimdVector tdb,
Simd< bool > &  rz_outside 
) const

Get complete B-field info

Parameters
[in]xVector of {r,z} coordinates
[in]fld_zVector of z coordinates
[out]bvecVector of magnetic field components (r,z,phi)
[out]jacbVector of jacobians (3x3xSIMD_SIZE)
[out]psivecVector of psi coordinates
[out]gradpsiVector of psi gradients (r,z)
[out]tdbVector of perturbed magnetic field (r,z,phi) NOT functional yet
[out]rz_outsideVector of whether the particles are outside the equilibrium

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::follow_field ( const SimdVector2D x_org,
const Simd< double > &  phi_org,
const Simd< double > &  phi_dest,
SimdVector2D x_dest 
) const

For a vector of locations in (r,z,phi), follow the equilibrium magnetic field back to a phi midplane to do particle-grid operations. The field is followed with a Runge Kutta scheme, which can be order 1, 2, or 4. Splitting this up to avoid branching could improve vectorization.

Parameters
[in]x_orgVector of (r,z) coordinates
[in]phi_orgVector of phi coordinates
[in]phi_destVector of phi coordinates of the target midplanes
[out]x_destVector of the (r,z) coordinates when the particles are projected along the field onto the specified plane

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::geometry_r ( double  r) const

Returns the first input argument if using normal toroidal geometry, and the r coordinate of the magnetic axis if cylindrical

Parameters
[in]ris the actual r coordinate
Returns
double the r suitable for the geometry

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi ( const SimdVector2D x,
Simd< double > &  psi_out 
) const

Get psi

Parameters
[in]r_inVector of r coordinates
[in]z_inVector of z coordinates
[out]psi_outVector of psi coordinates

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::get_psi ( double  r,
double  z 
) const
template<class Device >
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi_unit_vec ( const SimdVector2D x,
Simd< double > &  cosa,
Simd< double > &  sina 
) const

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::I_deriv ( double  psi_in,
int  rgn3 
) const

Interpolate the phi component of the magnetic field or its derivative at a single point

Parameters
[in]psi_inPoloidal magnetic flux at that point
[in]rgn3Whether the location is in region 3
Returns
The phi component of the magnetic field or its derivative
template<class Device >
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::I_value ( double  psi_in,
int  rgn3 
) const

Interpolate the phi component of the magnetic field or its derivative at a single point

Parameters
[in]psi_inPoloidal magnetic flux at that point
[in]rgn3Whether the location is in region 3
Returns
The phi component of the magnetic field or its derivative
template<class Device>
template<class Device2 >
MagneticField<Device2> MagneticField< Device >::mirror ( ) const
inline

Member Data Documentation

template<class Device>
RZBounds MagneticField< Device >::bounds

Simulation boundary.

template<class Device>
double MagneticField< Device >::bp_sign

Whether poloidal field is reversed?

template<class Device>
double MagneticField< Device >::bt_sign

Whether toroidal field is reversed?

template<class Device>
Equilibrium MagneticField< Device >::equil

The object containing information about the magnetic equilibrium.

template<class Device>
int MagneticField< Device >::ff_order

Order of RK scheme used for field following. Can be 1, 2, or 4.

template<class Device>
int MagneticField< Device >::ff_step

Number of steps taken when projecting the particle location onto the midplane.

template<class Device>
CubInterp<Device> MagneticField< Device >::I_interp

The object for interpolating I (for deriving toroidal magnetic field)

template<class Device>
double MagneticField< Device >::inpsi

Boundary condition used in a few spots.

template<class Device>
double MagneticField< Device >::outpsi

Boundary condition used in a few spots.

template<class Device>
Bicub<Device> MagneticField< Device >::psi_bicub

The object for interpolating psi (magnetic flux surfaces)


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