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

#include <grid.hpp>

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

Public Member Functions

 Grid (NLReader::NamelistReader &nlr, const GridFiles &grid_files, const MagneticField< DeviceType > &magnetic_field)
 
 Grid ()
 
template<class Device2 >
Grid< Device2 > mirror () const
 
KOKKOS_INLINE_FUNCTION void t_coeff (const SimdVector2D &x, const Simd< int > &itr, SimdGridVec &p) const
 
KOKKOS_INLINE_FUNCTION void t_coeff_mod (const MagneticField< Device > &magnetic_field, const double r, const double z, const double psiin, const int itr, double(&p)[3]) const
 
KOKKOS_INLINE_FUNCTION void t_coeff_mod (const MagneticField< Device > &magnetic_field, const SimdVector2D &xy, const Simd< double > &psiin, const Simd< int > &itr, SimdGridVec &p) const
 
KOKKOS_INLINE_FUNCTION void search_tr2_no_precheck (const double r, const double z, int &itr, double(&p)[3]) const
 
KOKKOS_INLINE_FUNCTION void search_tr2 (const SimdVector2D &xy, Simd< int > &itr, SimdGridVec &pout) const
 
KOKKOS_INLINE_FUNCTION void search_tr_check_guess (const SimdVector2D &x, const Simd< int > &old_itr, Simd< int > &itr, SimdGridVec &p) const
 
KOKKOS_INLINE_FUNCTION void psi_search (double psi, double &wp, int &ip) const
 
KOKKOS_INLINE_FUNCTION void follow_field_to_midplane (const MagneticField< Device > &magnetic_field, const SimdVector2D &x, const Simd< double > &phi, SimdVector2D &xff) const
 
KOKKOS_INLINE_FUNCTION void get_grid_weights_ff (const MagneticField< Device > &magnetic_field, const SimdVector &v, const Simd< double > &psi_in, SimdVector2D &xff, SimdGridWeights< Order::One, PhiInterpType::Planes > &grid_wts) const
 
KOKKOS_INLINE_FUNCTION void get_grid_weights_no_ff (const MagneticField< Device > &magnetic_field, const SimdVector2D &x, const Simd< double > &psi_in, SimdGridWeights< Order::One, PhiInterpType::None > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void nearest_node (const SimdGridWeights< Order::One, PIT > &grid_wts, SimdGridWeights< Order::Zero, PIT > &grid_wts0) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, const Simd< double > &psi, SimdVector2D &xff, SimdGridWeights< Order::One, PIT > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, SimdVector2D &xff, SimdGridWeights< Order::One, PIT > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, const Simd< double > &psi_in, SimdGridWeights< Order::Zero, PIT > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, const Simd< double > &psi_in, SimdGridWeights< Order::One, PIT > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, SimdGridWeights< Order::Zero, PIT > &grid_wts) const
 
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void get_grid_weights (const MagneticField< Device > &magnetic_field, const SimdVector &v, SimdGridWeights< Order::One, PIT > &grid_wts) const
 
KOKKOS_INLINE_FUNCTION void get_wall_index (const Simd< bool > &just_left_the_grid, const SimdVector2D &x, const SimdGridWeights< Order::One, PIT_GLOBAL > &grid_wts, Simd< int > &widx) const
 
KOKKOS_INLINE_FUNCTION void wedge_modulo_phi (Simd< double > &phi_mod) const
 
KOKKOS_INLINE_FUNCTION double wedge_modulo_phi (double phi) const
 
KOKKOS_INLINE_FUNCTION bool node_is_inside_psi_range (const MagneticField< Device > &magnetic_field, const int node) const
 
KOKKOS_INLINE_FUNCTION int get_plane_index (double phi) const
 
KOKKOS_INLINE_FUNCTION int get_node_index (int triangle_index, int tri_vertex_index) const
 
KOKKOS_INLINE_FUNCTION double node_area_to_volume (const MagneticField< Device > &magnetic_field, double area, int node_index) const
 
KOKKOS_INLINE_FUNCTION void get_triangle_area_and_volume (const MagneticField< Device > &magnetic_field, int i, double &area, double &volume) const
 
KOKKOS_INLINE_FUNCTION double get_r_center_of_mass (int itr) const
 
KOKKOS_INLINE_FUNCTION bool node_is_in_region_1_or_2_no_wall (const int inode) const
 
KOKKOS_INLINE_FUNCTION bool node_is_in_private_region_no_wall (const int inode) const
 
KOKKOS_INLINE_FUNCTION bool node_is_in_included_region (const int inode, const bool exclude_private_region) const
 
KOKKOS_INLINE_FUNCTION bool node_is_in_region_1_or_2 (const MagneticField< Device > &magnetic_field, const int inode) const
 
KOKKOS_INLINE_FUNCTION int uses_rz_basis (const int inode) const
 
KOKKOS_INLINE_FUNCTION double get_r (const int inode) const
 
KOKKOS_INLINE_FUNCTION double get_dist2_from_node (const int inode, const double r, const double z) const
 
KOKKOS_INLINE_FUNCTION double get_dist_from_node (const int inode, const double r, const double z) const
 
KOKKOS_INLINE_FUNCTION void get_rz_coordinates (const int inode, double &r, double &z) const
 
KOKKOS_INLINE_FUNCTION void get_rz_coordinates (const Simd< int > &grid_inds, SimdVector2D &x) const
 
KOKKOS_INLINE_FUNCTION void get_rz_coordinates (const int itr, const double(&p)[3], double &r, double &z) const
 
KOKKOS_INLINE_FUNCTION void get_rz_coordinates (const Simd< int > &itr, const SimdGridVec &p, SimdVector2D &x) const
 
KOKKOS_INLINE_FUNCTION int get_nearest_node (const Simd< int > &itr, const SimdGridVec &p, int i_simd) const
 
KOKKOS_INLINE_FUNCTION int get_nearest_node (const int itr, const double(&p)[3]) const
 
KOKKOS_INLINE_FUNCTION double get_nearest_midplane (double phi) const
 
KOKKOS_INLINE_FUNCTION void set_gradient_mat_triangle (const View< int *, CLayout, DeviceType > &num_t_node, const View< int **, CLayout, DeviceType > &tr_node, const View< double *, CLayout, DeviceType > &tr_area, const View< double ***, CLayout, DeviceType > &unit_vecs, const int i, const GradientMatrices< Device > &gradient_matrices) const
 
KOKKOS_INLINE_FUNCTION int get_characteristic_length (const View< int *, CLayout, DeviceType > &num_t_node, const View< int **, CLayout, DeviceType > &tr_node, const int i, double &dist_triangle, double &dist_psi, double &dist_theta) const
 
KOKKOS_INLINE_FUNCTION void set_grad_matrix_from_psi_theta (bool is_psi_dir, int itr_pos, int itr_neg, double(&p_pos)[3], double(&p_neg)[3], double dl_pos, double dl_neg, int inode, const Matrix< Device > &matrix) const
 
void draw_ascii_grid (MagneticField< Device > magnetic_field, double rmin, double rmax, double dr, double zmin, double zmax, double dz)
 

Public Attributes

int ntriangle
 Number of grid triangles. More...
 
int nnode
 Number of grid nodes. More...
 
int nplanes
 Number of planes. More...
 
double wedge_angle
 The size of the wedge (the model is periodic in phi, a angle of e.g. pi means half the tokamak is modeled) More...
 
double delta_phi
 Distance between planes. More...
 
double inv_delta_phi
 1/delta_phi More...
 
Kokkos::View< double
*, Kokkos::LayoutRight, Device > 
psi
 An array of psi coordinates. More...
 
Kokkos::View< double
*, Kokkos::LayoutRight, Device > 
bfield
 Magnetic field magnitude at nodes. More...
 
int npsi_surf2
 
double minval_psi_surf2
 
double maxval_psi_surf2
 
Kokkos::View< double
*, Kokkos::LayoutRight, Device > 
psi_surf2
 
UniformRange psi00
 
int nwall
 
VolumesAndAreas volumes_and_areas
 
Projection< HostTypehalf_plane_ff
 
Projection< HostTypeone_plane_ff
 

Private Attributes

Kokkos::View< int
*, Kokkos::LayoutRight, Device > 
basis
 0 for Psi-theta, 1 for R-Z ; should be enum (and bool/char?) More...
 
Kokkos::View< Vertex
*, Kokkos::LayoutRight, Device > 
nodes
 
Kokkos::View< VertexMap
*, Kokkos::LayoutRight, Device > 
mapping
 
GuessTable< Device > guess
 
GuessList1D< Device > psi_guess
 
Kokkos::View< int
*, Kokkos::LayoutRight, Device > 
rgn
 
Kokkos::View< RZPair
*, Kokkos::LayoutRight, Device > 
gx
 
View< int **, CLayout, Device > adj
 The indices of the three adjacent triangles for each triangle. More...
 
Kokkos::View< int
*, Kokkos::LayoutRight, Device > 
wall_nodes
 
Kokkos::View< int
*, Kokkos::LayoutRight, Device > 
node_to_wall
 
double eps_flux_surface
 The minimum difference in psi that counts as being on two different flux surfaces. More...
 

Friends

template<class Device2 >
class Grid
 

Constructor & Destructor Documentation

template<class Device >
Grid< Device >::Grid ( NLReader::NamelistReader nlr,
const GridFiles grid_files,
const MagneticField< DeviceType > &  magnetic_field 
)

< Simulate a wedge of 2pi/sml_wedge_n of the full torus

< Size of the hash table used to narrow down the number of triangles in which to search a particle

< Calculate node volume using Monte-Carlo method (true), or analytic method (false)

< If \(|\nabla\psi|\) is smaller than this threshold, \(\nabla\psi\) is computed

< Size of the uniform 1D psi-grid

Here is the call graph for this function:

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

Member Function Documentation

template<class Device>
void Grid< Device >::draw_ascii_grid ( MagneticField< Device >  magnetic_field,
double  rmin,
double  rmax,
double  dr,
double  zmin,
double  zmax,
double  dz 
)
inline
template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::follow_field_to_midplane ( const MagneticField< Device > &  magnetic_field,
const SimdVector2D x,
const Simd< double > &  phi,
SimdVector2D xff 
) const

Follow the magnetic field to the midplane

Parameters
[in]magnetic_fieldThe magnetic field object
[in]xVector of (r,z) coordinates
[in]phiVector of phi coordinates
[out]xffVector of (r,z) coordinates mapped along B-field to midplane

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::get_characteristic_length ( const View< int *, CLayout, DeviceType > &  num_t_node,
const View< int **, CLayout, DeviceType > &  tr_node,
const int  i,
double &  dist_triangle,
double &  dist_psi,
double &  dist_theta 
) const

Determines characteristic distance in grad(psi) and grad(theta) direction

Parameters
[in]num_t_node
[in]tr_node
[in]inode index
[out]dist_triangle
[out]dist_psi
[out]dist_theta
Returns
int error code: -3 is fatal, -2 means all vertices are on the same surface, -1 means?

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::get_dist2_from_node ( const int  inode,
const double  r,
const double  z 
) const

Gets the square of the distance from a grid vertex of the input (r,z) coordinates. The square of the distance is sufficient in some algorithms, and avoids a costly square root. If the distance itself is needed, use get_dist_from_node

Parameters
[in]inodegrid node
[in]rr coordinate of input point
[in]zz coordinate of input point
Returns
double the square of the distance

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::get_dist_from_node ( const int  inode,
const double  r,
const double  z 
) const

Gets the distance from a grid vertex of the input (r,z) coordinates.

Parameters
[in]inodegrid node
[in]rr coordinate of input point
[in]zz coordinate of input point
Returns
double the distance

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
const Simd< double > &  psi_in,
SimdVector2D xff,
SimdGridWeights< Order::One, PIT > &  grid_wts 
) const

Get the grid triangle and vertex weights for a vector of locations. Name is historical and should change

Parameters
[in]magnetic_fieldThe magnetic field object (currently global so technically unnecessary)
[in]vVector of (r,z, phi) coordinates
[in]phiVector of phi coordinates
[in]psi_inVector of psi coordinates
[out]xffVector of (r,z) coordinates mapped along B-field to midplane
[out]itrVector of grid triangles
[out]pVector of triangle vertex weights

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
SimdVector2D xff,
SimdGridWeights< Order::One, PIT > &  grid_wts 
) const
template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
const Simd< double > &  psi_in,
SimdGridWeights< Order::Zero, PIT > &  grid_wts 
) const

Here is the call graph for this function:

template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
const Simd< double > &  psi,
SimdGridWeights< Order::One, PIT > &  grid_wts 
) const

Get the grid triangle and vertex weights for a vector of locations. Name is historical and should change

Parameters
[in]magnetic_fieldThe magnetic field object
[in]vVector of (r,z,phi) coordinates
[in]psiVector of psi coordinates
[out]itrVector of grid triangles
[out]pVector of triangle vertex weights

Here is the call graph for this function:

template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
SimdGridWeights< Order::Zero, PIT > &  grid_wts 
) const

Here is the call graph for this function:

template<class Device>
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
SimdGridWeights< Order::One, PIT > &  grid_wts 
) const

Get the grid triangle and vertex weights for a vector of locations. Name is historical and should change

Parameters
[in]magnetic_fieldThe magnetic field object
[in]vVector of (r,z,phi) coordinates
[out]itrVector of grid triangles
[out]pVector of triangle vertex weights

Here is the call graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights_ff ( const MagneticField< Device > &  magnetic_field,
const SimdVector v,
const Simd< double > &  psi_in,
SimdVector2D xff,
SimdGridWeights< Order::One, PhiInterpType::Planes > &  grid_wts 
) const

Get the grid triangle and vertex weights for a vector of locations. Name is historical and should change

Parameters
[in]magnetic_fieldThe magnetic field object
[in]xVector of (r,z) coordinates
[in]phiVector of phi coordinates
[in]psi_inVector of psi coordinates
[out]xffVector of (r,z) coordinates mapped along B-field to midplane
[out]itrVector of grid triangles
[out]pVector of triangle vertex weights

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_grid_weights_no_ff ( const MagneticField< Device > &  magnetic_field,
const SimdVector2D x,
const Simd< double > &  psi_in,
SimdGridWeights< Order::One, PhiInterpType::None > &  grid_wts 
) const

Get the grid triangle and vertex weights for a vector of locations. Name is historical and should change

Parameters
[in]magnetic_fieldThe magnetic field object
[in]xVector of (r,z) coordinates
[in]phiVector of phi coordinates
[in]psi_inVector of psi coordinates
[out]xffVector of (r,z) coordinates mapped along B-field to midplane
[out]itrVector of grid triangles
[out]pVector of triangle vertex weights

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::get_nearest_midplane ( double  phi) const

Returns nearest midplane for a given phi

Parameters
[in]phiinput phi coordinate
template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::get_nearest_node ( const Simd< int > &  itr,
const SimdGridVec p,
int  i_simd 
) const

Returns nearest node given a triangle and barycentric coordinates

Parameters
[in]itrtriangle index (1-indexed)
[in]pbarycentric weights
[in]i_simdsimd index
Returns
nearest node (1-indexed)
template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::get_nearest_node ( const int  itr,
const double(&)  p[3] 
) const

Returns nearest node given a triangle and barycentric coordinates

Parameters
[in]itrtriangle index (1-indexed)
[in]pbarycentric weights
Returns
nearest node (1-indexed)
template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::get_node_index ( int  triangle_index,
int  tri_vertex_index 
) const

Returns the index of the node (0-indexed), for a given triangle and triangle vertex index

Parameters
[in]triangle_indexthe index of the triangle (1-indexed)
[in]tri_vertex_indexthe index of the requested vertex on the triangle
Returns
int the index of the node (0-indexed)

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::get_plane_index ( double  phi) const

Returns the index of the plane, rounded down.

Parameters
[in]phithe phi value
Returns
int the index of the plane

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::get_r ( const int  inode) const

Accesses the r coordinate of the grid node

Parameters
[in]inodegrid node to check
Returns
double the r coordinate

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::get_r_center_of_mass ( int  itr) const

Returns r coordinate of center of mass of a triangle

Parameters
[in]iis the triangle index (0-indexed)
Returns
double r_center_of_mass
template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_rz_coordinates ( const int  inode,
double &  r,
double &  z 
) const

Returns rz coordinates of a single node

Parameters
[in]inodenode index
[out]rr coordinate
[out]zz coordinate

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_rz_coordinates ( const Simd< int > &  grid_inds,
SimdVector2D x 
) const

Returns vector of rz coordinates of the input node indices

Parameters
[in]grid_indsVector of node indices
[out]xVector of rz coordinates
template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_rz_coordinates ( const int  itr,
const double(&)  p[3],
double &  r,
double &  z 
) const

Returns rz coordinates given a triangle index and barycentric weights

Parameters
[in]itrtriangle index (1-indexed)
[in]pbarycentric weights
[out]rr coordinate of input point
[out]zz coordinate of input point

Here is the call graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_rz_coordinates ( const Simd< int > &  itr,
const SimdGridVec p,
SimdVector2D x 
) const

Returns rz coordinates given a triangle index and barycentric weights

Parameters
[in]itrtriangle index (1-indexed)
[in]pbarycentric weights
[out]xr,z coordinates of input point

Here is the call graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_triangle_area_and_volume ( const MagneticField< Device > &  magnetic_field,
int  i,
double &  area,
double &  volume 
) const

Returns analytical calculation of triangle area and projected volume

Parameters
[in]magnetic_fieldis the magnetic field
[in]iis the triangle index (0-indexed)
[out]areais the triangle area
[out]volumeis the triangle projected volume
Returns
void

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::get_wall_index ( const Simd< bool > &  just_left_the_grid,
const SimdVector2D x,
const SimdGridWeights< Order::One, PIT_GLOBAL > &  grid_wts,
Simd< int > &  widx 
) const

Find nearest wall point for a vector of locations.

Parameters
[in]just_left_the_gridVector of whether particles have left the grid
[in]xVector of (r,z) coordinates
[in]itrVector of grid triangles the particles were last in
[in]pVector of the vertex weighting in those triangles
[out]widxVector of wall point indices

Here is the caller graph for this function:

template<class Device>
template<class Device2 >
Grid<Device2> Grid< Device >::mirror ( ) const
inline
template<class Device >
template<PhiInterpType PIT>
KOKKOS_INLINE_FUNCTION void Grid< Device >::nearest_node ( const SimdGridWeights< Order::One, PIT > &  grid_wts,
SimdGridWeights< Order::Zero, PIT > &  grid_wts0 
) const

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION double Grid< Device >::node_area_to_volume ( const MagneticField< Device > &  magnetic_field,
double  area,
int  node_index 
) const

Returns analytical calculation of node volume by projecting the node area into the toroidal direction Assumes the center of mass is the node

Parameters
[in]magnetic_fieldis the magnetic field
[in]areais the area of the node
[in]node_indexis the node index
Returns
double the area owned by the node

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION bool Grid< Device >::node_is_in_included_region ( const int  inode,
const bool  exclude_private_region 
) const

Checks if a grid node lies within an included region. This excludes wall nodes and the private region if specified.

Parameters
[in]inodegrid node to check
[in]exclude_private_regionwhether to private region is excluded
Returns
bool true if in included region

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION bool Grid< Device >::node_is_in_private_region_no_wall ( const int  inode) const

Checks if a grid node lies within region 3. Excludes wall since rgn is a different value there

Parameters
[in]inodegrid node to check
Returns
bool true if in region 3

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION bool Grid< Device >::node_is_in_region_1_or_2 ( const MagneticField< Device > &  magnetic_field,
const int  inode 
) const

Checks if a grid node lies within region 1 or 2. Note this is subtly different from rgn, because rgn also has a wall option

Parameters
[in]magnetic_fieldthe magnetic field object
[in]inodegrid node to check
Returns
bool true if in region 1 or 2

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION bool Grid< Device >::node_is_in_region_1_or_2_no_wall ( const int  inode) const

Checks if a grid node lies within region 1 or 2. Excludes wall since rgn is a different value there

Parameters
[in]inodegrid node to check
Returns
bool true if in region 1 or 2

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION bool Grid< Device >::node_is_inside_psi_range ( const MagneticField< Device > &  magnetic_field,
const int  node 
) const

Checks if a grid node lies within the simulation domain (radially).

Parameters
[in]magnetic_fieldthe magnetic field object (currently global so technically unnecessary)
[in]nodegrid node to check (0-indexed)
Returns
bool true if inside, false if outside

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::psi_search ( double  psi,
double &  wp,
int &  ip 
) const

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::search_tr2 ( const SimdVector2D xy,
Simd< int > &  itr,
SimdGridVec pout 
) const

Determine triangle location and triangle vertex weights of vector of locations, first checking the particles have left their previous triangle, and if so doing a search

Parameters
[in]xyVector of (r,z) coordinates
[out]itrVector of triangle locations
[out]poutVector of triangle vertex weights

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::search_tr2_no_precheck ( const double  r,
const double  z,
int &  itr,
double(&)  p[3] 
) const

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::search_tr_check_guess ( const SimdVector2D x,
const Simd< int > &  old_itr,
Simd< int > &  itr,
SimdGridVec p 
) const

A check to see if vector of particles has stayed in the same triangle since last time. If so just recalculate p weights

Parameters
[in]xVector of (r,z) coordinates
[in]old_itrVector of old triangle location
[out]itrVector of new triangle location (same as old if unchanged, -1 if search is needed)
[out]pVector of p weights if in same triangle

Here is the call graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::set_grad_matrix_from_psi_theta ( bool  is_psi_dir,
int  itr_pos,
int  itr_neg,
double(&)  p_pos[3],
double(&)  p_neg[3],
double  dl_pos,
double  dl_neg,
int  inode,
const Matrix< Device > &  matrix 
) const

Sets gradient matrix values based on psi-theta triangle weights

Parameters
[in]is_psi_dirwhether the gradient matrix is psi or theta direction (just for the error message)
[in]itr_posthe triangle index in the positive direction
[in]itr_negthe triangle index in the negative direction
[in]p_posthe triangle weights in the positive direction
[in]p_negthe triangle weights in the negative direction
[in]dl_posthe distance in the positive direction
[in]dl_negthe distance in the negative direction
[in]inodethe node index
[out]matrixthe gradient matrix being updates

Here is the call graph for this function:

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::set_gradient_mat_triangle ( const View< int *, CLayout, DeviceType > &  num_t_node,
const View< int **, CLayout, DeviceType > &  tr_node,
const View< double *, CLayout, DeviceType > &  tr_area,
const View< double ***, CLayout, DeviceType > &  unit_vecs,
const int  i,
const GradientMatrices< Device > &  gradient_matrices 
) const

Sets gradient matrix values for a given vertex

Parameters
[in]num_t_node
[in]tr_node
[in]tr_area
[in]unit_vecs
[in]inode index
[out]gradient_matricesGradient matrices

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::t_coeff ( const SimdVector2D x,
const Simd< int > &  itr,
SimdGridVec p 
) const

Calculate a vector of barycentric coordinates p from real-space coordinates x in triangle itr.

Parameters
[in]xVector of (r,z) coordinates
[in]itrVector of triangles
[out]pVector of the weights of the triangle vertices

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::t_coeff_mod ( const MagneticField< Device > &  magnetic_field,
const double  r,
const double  z,
const double  psiin,
const int  itr,
double(&)  p[3] 
) const

Here is the caller graph for this function:

template<class Device>
KOKKOS_INLINE_FUNCTION void Grid< Device >::t_coeff_mod ( const MagneticField< Device > &  magnetic_field,
const SimdVector2D xy,
const Simd< double > &  psiin,
const Simd< int > &  itr,
SimdGridVec p 
) const
template<class Device >
KOKKOS_INLINE_FUNCTION int Grid< Device >::uses_rz_basis ( const int  inode) const

Whether the field uses an r-z or a psi-theta basis.

Parameters
[in]inodegrid node to check
Returns
int 1 if rz basis, 0 if psi-theta basis

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION void Grid< Device >::wedge_modulo_phi ( Simd< double > &  phi_mod) const

Modulo phi in place: Keep phi between 0 and wedge_angle. Best algorithm as long as the particles haven't left the wedge multiple times over

Parameters
[in,out]phi_modVector of phi coordinates

Here is the caller graph for this function:

template<class Device >
KOKKOS_INLINE_FUNCTION double Grid< Device >::wedge_modulo_phi ( double  phi) const

Modulo phi in place: Keep phi between 0 and wedge_angle. Best algorithm as long as the particles haven't left the wedge multiple times over

Parameters
[in,out]phi_modVector of phi coordinates

Friends And Related Function Documentation

template<class Device>
template<class Device2 >
friend class Grid
friend

Member Data Documentation

template<class Device>
View<int**,CLayout,Device> Grid< Device >::adj
private

The indices of the three adjacent triangles for each triangle.

template<class Device>
Kokkos::View<int*, Kokkos::LayoutRight,Device> Grid< Device >::basis
private

0 for Psi-theta, 1 for R-Z ; should be enum (and bool/char?)

template<class Device>
Kokkos::View<double*,Kokkos::LayoutRight,Device> Grid< Device >::bfield

Magnetic field magnitude at nodes.

template<class Device>
double Grid< Device >::delta_phi

Distance between planes.

template<class Device>
double Grid< Device >::eps_flux_surface
private

The minimum difference in psi that counts as being on two different flux surfaces.

template<class Device>
GuessTable<Device> Grid< Device >::guess
private
template<class Device>
Kokkos::View<RZPair*,Kokkos::LayoutRight,Device> Grid< Device >::gx
private
template<class Device>
Projection<HostType> Grid< Device >::half_plane_ff
template<class Device>
double Grid< Device >::inv_delta_phi

1/delta_phi

template<class Device>
Kokkos::View<VertexMap*,Kokkos::LayoutRight,Device> Grid< Device >::mapping
private
template<class Device>
double Grid< Device >::maxval_psi_surf2
template<class Device>
double Grid< Device >::minval_psi_surf2
template<class Device>
int Grid< Device >::nnode

Number of grid nodes.

template<class Device>
Kokkos::View<int*,Kokkos::LayoutRight,Device> Grid< Device >::node_to_wall
private
template<class Device>
Kokkos::View<Vertex*,Kokkos::LayoutRight,Device> Grid< Device >::nodes
private
template<class Device>
int Grid< Device >::nplanes

Number of planes.

template<class Device>
int Grid< Device >::npsi_surf2
template<class Device>
int Grid< Device >::ntriangle

Number of grid triangles.

template<class Device>
int Grid< Device >::nwall
template<class Device>
Projection<HostType> Grid< Device >::one_plane_ff
template<class Device>
Kokkos::View<double*,Kokkos::LayoutRight,Device> Grid< Device >::psi

An array of psi coordinates.

template<class Device>
UniformRange Grid< Device >::psi00
template<class Device>
GuessList1D<Device> Grid< Device >::psi_guess
private
template<class Device>
Kokkos::View<double*,Kokkos::LayoutRight,Device> Grid< Device >::psi_surf2
template<class Device>
Kokkos::View<int*,Kokkos::LayoutRight,Device> Grid< Device >::rgn
private
template<class Device>
VolumesAndAreas Grid< Device >::volumes_and_areas
template<class Device>
Kokkos::View<int*,Kokkos::LayoutRight,Device> Grid< Device >::wall_nodes
private
template<class Device>
double Grid< Device >::wedge_angle

The size of the wedge (the model is periodic in phi, a angle of e.g. pi means half the tokamak is modeled)


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