XGCa
|
#include <magnetic_field.hpp>
Public Member Functions | |
MagneticField (NLReader::NamelistReader &nlr, const MagneticEquilFiles::Ptr &equil_files, const View< Equil::XPoint *, HostType > &xpts, const RZPair &axis_in, const RZPair &psi_ref_point) | |
void | write () const |
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 double | psi_norm () const |
KOKKOS_INLINE_FUNCTION void | get_psi (const SimdVector &v, Simd< double > &psi_out) const |
KOKKOS_INLINE_FUNCTION void | get_psi (const SimdVector2D &x, const Simd< double > &phi, Simd< double > &psi_out) const |
KOKKOS_INLINE_FUNCTION double | get_psi (double r, double z, double phi) const |
KOKKOS_INLINE_FUNCTION void | get_psi_and_derivs (double r, double z, double phi, double &psi, double &dpsidr, double &dpsidz, double &dpsidphi) const |
KOKKOS_INLINE_FUNCTION void | get_psi_unit_vec (const SimdVector2D &x, double phi, Simd< double > &cosa, Simd< double > &sina) const |
KOKKOS_INLINE_FUNCTION double | geometry_r (double r) const |
KOKKOS_INLINE_FUNCTION bool | is_in_region_1_or_2 (double r, double z, double psi) const |
KOKKOS_INLINE_FUNCTION bool | is_in_region_1 (double r, double z, double psi) const |
KOKKOS_INLINE_FUNCTION void | check_boundaries (const SimdVector2D &x, Simd< bool > &rz_outside) const |
KOKKOS_INLINE_FUNCTION void | get_theta (const SimdVector2D &x, Simd< double > &theta) const |
KOKKOS_INLINE_FUNCTION void | field (const SimdVector &v, 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<> | |
void | write () const |
Public Attributes | |
double | bt_sign |
Whether toroidal field is reversed? More... | |
double | bp_sign |
Whether poloidal field is reversed? 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... | |
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 |
KOKKOS_INLINE_FUNCTION double | I_value (double psi_in, int rgn3) const |
KOKKOS_INLINE_FUNCTION double | I_deriv (double psi_in, int rgn3) const |
template<> | |
MagneticField (NLReader::NamelistReader &nlr, const MagneticEquilFiles::Ptr &equil_files_ptr, const View< Equil::XPoint *, HostType > &xpts, const RZPair &axis_in, const RZPair &psi_ref_point) | |
template<> | |
MagneticField (PsiOption psi_opt, double safety_factor_coeff, int eq_mr_in, int eq_mz_in, int eq_mpsi_in) | |
Private Attributes | |
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... | |
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... | |
Friends | |
template<class Device2 > | |
class | MagneticField |
MagneticField< Device >::MagneticField | ( | NLReader::NamelistReader & | nlr, |
const MagneticEquilFiles::Ptr & | equil_files, | ||
const View< Equil::XPoint *, HostType > & | xpts, | ||
const RZPair & | axis_in, | ||
const RZPair & | psi_ref_point | ||
) |
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 |
|
private |
Constructor for magnetic field
|
private |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::bmag_interpol | ( | const SimdVector & | v, |
Simd< double > & | bmag | ||
) | const |
For vector of coordinates, get the B-field magnitude
[in] | v | Vector of {r,z,phi} coordinates |
[out] | bmag | Vector of magnetic field magnitudes |
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
[in] | x | Vector of {r,z} coordinates |
[out] | bmag | Vector of magnetic field magnitudes |
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
[in] | r | r coordinate |
[in] | z | z coordinate |
[in] | phi | phi coordinate |
[out] | br | r component of B |
[out] | bz | z component of B |
[out] | bphi | phi component of B |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::check_boundaries | ( | const SimdVector2D & | x, |
Simd< bool > & | rz_outside | ||
) | const |
|
private |
Get the derivatives of the magnetic field to see where the particle should move to as it traces the field
[in] | x | (r,z) coordinates |
[in] | phi | phi coordinate |
[out] | dx | change in (r,z) |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::field | ( | const SimdVector & | v, |
SimdVector & | bvec, | ||
SimdVector(&) | jacb[3], | ||
Simd< double > & | psivec, | ||
SimdVector2D & | gradpsi, | ||
SimdVector & | tdb, | ||
Simd< bool > & | rz_outside | ||
) | const |
Get complete B-field info
[in] | v | Vector of {r,z, phi} coordinates |
[out] | bvec | Vector of magnetic field components (r,z,phi) |
[out] | jacb | Vector of jacobians (3x3xSIMD_SIZE) |
[out] | psivec | Vector of psi coordinates |
[out] | gradpsi | Vector of psi gradients (r,z) |
[out] | tdb | Vector of perturbed magnetic field (r,z,phi) NOT functional yet |
[out] | rz_outside | Vector of whether the particles are outside the equilibrium |
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.
[in] | x_org | Vector of (r,z) coordinates |
[in] | phi_org | Vector of phi coordinates |
[in] | phi_dest | Vector of phi coordinates of the target midplanes |
[out] | x_dest | Vector of the (r,z) coordinates when the particles are projected along the field onto the specified plane |
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
[in] | r | is the actual r coordinate |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi | ( | const SimdVector & | v, |
Simd< double > & | psi_out | ||
) | const |
Get psi
[in] | r_in | Vector of r coordinates |
[in] | z_in | Vector of z coordinates |
[out] | psi_out | Vector of psi coordinates |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi | ( | const SimdVector2D & | x, |
const Simd< double > & | phi, | ||
Simd< double > & | psi_out | ||
) | const |
Get psi
[in] | r_in | Vector of r coordinates |
[in] | z_in | Vector of z coordinates |
[out] | psi_out | Vector of psi coordinates |
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::get_psi | ( | double | r, |
double | z, | ||
double | phi | ||
) | const |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi_and_derivs | ( | double | r, |
double | z, | ||
double | phi, | ||
double & | psi, | ||
double & | dpsidr, | ||
double & | dpsidz, | ||
double & | dpsidphi | ||
) | const |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_psi_unit_vec | ( | const SimdVector2D & | x, |
double | phi, | ||
Simd< double > & | cosa, | ||
Simd< double > & | sina | ||
) | const |
KOKKOS_INLINE_FUNCTION void MagneticField< Device >::get_theta | ( | const SimdVector2D & | x, |
Simd< double > & | theta | ||
) | const |
|
private |
Interpolate the phi component of the magnetic field or its derivative at a single point
[in] | psi_in | Poloidal magnetic flux at that point |
[in] | rgn3 | Whether the location is in region 3 |
|
private |
Interpolate the phi component of the magnetic field or its derivative at a single point
[in] | psi_in | Poloidal magnetic flux at that point |
[in] | rgn3 | Whether the location is in region 3 |
KOKKOS_INLINE_FUNCTION bool MagneticField< Device >::is_in_region_1 | ( | double | r, |
double | z, | ||
double | psi | ||
) | const |
KOKKOS_INLINE_FUNCTION bool MagneticField< Device >::is_in_region_1_or_2 | ( | double | r, |
double | z, | ||
double | psi | ||
) | const |
|
inline |
KOKKOS_INLINE_FUNCTION double MagneticField< Device >::psi_norm | ( | ) | const |
void MagneticField< HostType >::write | ( | ) | const |
void MagneticField< Device >::write | ( | ) | const |
RZBounds MagneticField< Device >::bounds |
Simulation boundary.
double MagneticField< Device >::bp_sign |
Whether poloidal field is reversed?
double MagneticField< Device >::bt_sign |
Whether toroidal field is reversed?
Equilibrium MagneticField< Device >::equil |
The object containing information about the magnetic equilibrium.
|
private |
Order of RK scheme used for field following. Can be 1, 2, or 4.
|
private |
Number of steps taken when projecting the particle location onto the midplane.
|
private |
The object for interpolating I (for deriving toroidal magnetic field)
double MagneticField< Device >::inpsi |
Boundary condition used in a few spots.
double MagneticField< Device >::outpsi |
Boundary condition used in a few spots.
|
private |
The object for interpolating psi (magnetic flux surfaces)