XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | Private Member Functions | 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 DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field)
 
 Grid ()
 
template<class Device2 >
Grid< Device2 > mirror () const
 
void write_to_file () const
 
KOKKOS_INLINE_FUNCTION void psi_search (double psi, double &wp, int &ip) 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 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 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 RZPair get_wall_rz (int i_wall) const
 
template<>
void write_to_file () const
 

Public Attributes

bool axisymmetric
 
Plane< Device > midplane
 
Plane< Device > lplane
 
Plane< Device > rplane
 
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...
 
int nwall
 
int npsi_surf2
 
double minval_psi_surf2
 
double maxval_psi_surf2
 
UniformRange psi00
 
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...
 
Kokkos::View< double
*, Kokkos::LayoutRight, Device > 
psi_surf2
 
Kokkos::View< int
*, Kokkos::LayoutRight, Device > 
wall_nodes
 
VolumesAndAreas volumes_and_areas
 
Projection< HostTypehalf_plane_ff
 
Projection< HostTypeff_rplane_to_neighbors
 
Projection< HostTypeff_lplane_to_neighbors
 
Projection< HostTypeff_to_midplane
 

Private Member Functions

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

Friends

template<class Device2 >
class Grid
 

Constructor & Destructor Documentation

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

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

Here is the call graph for this function:

Member Function Documentation

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 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
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
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
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
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 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_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 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 >
KOKKOS_INLINE_FUNCTION RZPair Grid< Device >::get_wall_rz ( int  i_wall) const

Returns the (r,z) coordinates of a wall vertex given the index in wall_nodes

Parameters
[in]i_wallis the wall index requested
Returns
RZPair the (r,z) coordinates of the wall vertex

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 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 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

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 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 caller graph for this function:

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
template<class Device>
void Grid< Device >::write_to_file ( ) const
template<>
void Grid< HostType >::write_to_file ( ) const

Here is the call graph for this function:

Friends And Related Function Documentation

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

Member Data Documentation

template<class Device>
bool Grid< Device >::axisymmetric
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>
Projection<HostType> Grid< Device >::ff_lplane_to_neighbors
template<class Device>
Projection<HostType> Grid< Device >::ff_rplane_to_neighbors
template<class Device>
Projection<HostType> Grid< Device >::ff_to_midplane
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>
Plane<Device> Grid< Device >::lplane
template<class Device>
double Grid< Device >::maxval_psi_surf2
template<class Device>
Plane<Device> Grid< Device >::midplane
template<class Device>
double Grid< Device >::minval_psi_surf2
template<class Device>
int Grid< Device >::nnode

Number of grid nodes.

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>
Kokkos::View<double*,Kokkos::LayoutRight,Device> Grid< Device >::psi

An array of psi coordinates.

template<class Device>
UniformRange Grid< Device >::psi00
template<class Device>
Kokkos::View<double*,Kokkos::LayoutRight,Device> Grid< Device >::psi_surf2
template<class Device>
Plane<Device> Grid< Device >::rplane
template<class Device>
VolumesAndAreas Grid< Device >::volumes_and_areas
template<class Device>
Kokkos::View<int*,Kokkos::LayoutRight,Device> Grid< Device >::wall_nodes
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: