Cardinal
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
NekRSMesh Class Reference

#include <NekRSMesh.h>

Inheritance diagram for NekRSMesh:
[legend]

Public Member Functions

 NekRSMesh (const InputParameters &parameters)
 
 NekRSMesh (const NekRSMesh &)=default
 
NekRSMeshoperator= (const NekRSMesh &other_mesh)=delete
 
virtual std::unique_ptr< MooseMesh > safeClone () const override
 
void saveInitialVolMesh ()
 Save the initial volumetric mesh for volume mirror-based moving mesh problems. More...
 
void initializePreviousDisplacements ()
 Initialize previous displacement values to zero for boundary mirror-based moving mesh problems. More...
 
int nekPolynomialOrder () const
 
int nBuildPerVolumeElem () const
 
int nBuildPerSurfaceElem () const
 
const std::vector< double > & nek_initial_x () const
 
const std::vector< double > & nek_initial_y () const
 
const std::vector< double > & nek_initial_z () const
 
std::vector< double > & prev_disp_x ()
 
std::vector< double > & prev_disp_y ()
 
std::vector< double > & prev_disp_z ()
 
const NekBoundaryCouplingboundaryCoupling () const
 
const NekVolumeCouplingvolumeCoupling () const
 
virtual void addElems ()
 Add all the elements in the mesh to the MOOSE data structures. More...
 
const order::NekOrderEnumorder () const
 
int numQuadraturePoints1D () const
 
int nekNumQuadraturePoints1D () const
 
const int & numElems () const
 Get the number of NekRS elements we rebuild in the MOOSE mesh. More...
 
const int & numVerticesPerElem () const
 Get the number of vertices per element in MOOSE's representation of nekRS's mesh. More...
 
int nodeIndex (const int gll_index) const
 Get the libMesh node index from nekRS's GLL index ordering. More...
 
const int & numSurfaceElems () const
 
const int & nekNumSurfaceElems () const
 
const int & numVerticesPerSurface () const
 
const int & numVolumeElems () const
 
const int & nekNumVolumeElems () const
 
const int & numVerticesPerVolume () const
 
const std::vector< int > * boundary () const
 
const bool & volume () const
 
Elem * boundaryElem () const
 Create a new element for a boundary mesh. More...
 
Elem * volumeElem () const
 Create a new element for a volume mesh. More...
 
virtual void buildMesh () override
 
virtual void buildDummyMesh ()
 
virtual void extractSurfaceMesh ()
 
virtual void extractVolumeMesh ()
 
int boundaryNodeIndex (const int gll_index) const
 
int volumeNodeIndex (const int gll_index) const
 
const Real & scaling () const
 
virtual void printMeshInfo () const
 Print diagnostic information related to the mesh. More...
 
int boundaryElemProcessorID (const int elem_id)
 
int volumeElemProcessorID (const int elem_id)
 
int facesOnBoundary (const int elem_id) const
 
bool exactMirror () const
 
std::vector< std::vector< int > > cornerIndices () const
 
int nMoosePerNek () const
 
void updateDisplacement (const int e, const double *src, const field::NekWriteEnum field)
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

void storeVolumeCoupling ()
 Store the rank-local element and rank ownership for volume coupling. More...
 
void storeBoundaryCoupling ()
 
int boundary_id (const int elem_id, const int face_id)
 
void faceVertices ()
 
void volumeVertices ()
 
void initializeMeshParams ()
 Initialize members for the mesh and determine the GLL-to-node mapping. More...
 

Protected Attributes

const bool & _volume
 Whether nekRS is coupled through volumes to MOOSE. More...
 
const std::vector< int > * _boundary
 Boundary ID(s) through which to couple Nek to MOOSE. More...
 
const order::NekOrderEnum _order
 Order of the surface interpolation between nekRS and MOOSE. More...
 
const bool & _exact
 
int _n_vertices_per_surface
 Number of vertices per surface. More...
 
int _n_vertices_per_volume
 Number of vertices per volume element. More...
 
const Real & _scaling
 Spatial scaling factor to apply to the mesh. More...
 
const unsigned int & _fluid_block_id
 Block ID for the fluid portion of the mesh mirror. More...
 
const unsigned int & _solid_block_id
 Block ID for the solid portion of the mesh mirror. More...
 
int _nek_polynomial_order
 Order of the nekRS solution. More...
 
int _n_surface_elems
 
int _n_build_per_surface_elem
 Number of MOOSE surface elements to build per NekRS surface element. More...
 
int _n_volume_elems
 
int _n_build_per_volume_elem
 Number of MOOSE volume elements to build per NekRS volume element. More...
 
int _n_elems
 
int _n_moose_per_nek
 
int(NekRSMesh::* _elem_processor_id )(const int elem_id)
 Function returning the processor id which should own each element. More...
 
int _n_vertices_per_elem
 Number of vertices per element, which depends on whether building a boundary/volume mesh. More...
 
std::vector< int > * _node_index
 Mapping of GLL nodes to libMesh nodes, which depends on whether building a boundary/volume mesh. More...
 
int _nek_n_surface_elems
 Total number of surface elements in the nekRS problem. More...
 
int _nek_n_volume_elems
 Total number of volume elements in the nekRS problem. More...
 
std::vector< int > _phase
 "Phase" for each element (fluid = 0, solid = 1) More...
 
std::vector< double > _x
 \(x\) coordinates of the current GLL points (which can move in time), for this rank More...
 
std::vector< double > _y
 \(y\) coordinates of the current GLL points (which can move in time), for this rank More...
 
std::vector< double > _z
 \(z\) coordinates of the current GLL points (which can move in time), for this rank More...
 
std::vector< double > _initial_x
 \(x\) coordinates of the initial GLL points, for this rank More...
 
std::vector< double > _initial_y
 \(y\) coordinates of the initial GLL points, for this rank More...
 
std::vector< double > _initial_z
 \(z\) coordinates of the initial GLL points, for this rank More...
 
std::vector< int > _bnd_node_index
 Mapping of boundary GLL indices to MooseMesh node indices. More...
 
std::vector< int > _vol_node_index
 Mapping of volume GLL indices to MooseMesh node indices. More...
 
std::vector< int > _side_index
 Mapping of side indices to libMesh side indices. More...
 
Elem *(NekRSMesh::* _new_elem )() const
 Function pointer to the type of new element to add. More...
 
NekBoundaryCoupling _boundary_coupling
 Data structure holding mapping information for boundary coupling. More...
 
NekVolumeCoupling _volume_coupling
 Data structure holding mapping information for volume coupling. More...
 
mesh_t * _nek_internal_mesh = nullptr
 Pointer to NekRS's internal mesh data structure. More...
 
std::vector< std::vector< int > > _corner_indices
 Corner indices for GLL points of mesh mirror elements. More...
 
std::vector< double > _prev_disp_x
 
std::vector< double > _prev_disp_y
 
std::vector< double > _prev_disp_z
 

Detailed Description

Representation of a nekRS surface mesh as a native MooseMesh. This is constructed by interpolating from the surface Gauss-Lobatto-Legendre points in nekRS to either a first-order (Quad4) or second-order (Quad9) mesh. This mesh is only constructed for a user-specified set of boundaries in the nekRS mesh with the 'boundary' parameter. Therefore, this class contains a mixture of information related to the nekRS mesh (that nekRS solves its equations on) versus the surface mesh constructed for data transfer with MOOSE (which is only used by nekRS for the purpose of transferring its solution). All information specific to the mesh nekRS actually uses for its solution are prefaced with either '_nek' or 'nek' to help with this distinction.

The nekRS mesh is currently implemented as a replicated mesh. On the nekRS side, an Allgather is used to get the surface geometry information on each nekRS process such that access from MOOSE can be performed on each process.

TODO: The extension to higher than a second-order representation requires some modifications to the formation of the mesh, as well as the interpolation matrices used in NekRSProblem, because for 3rd order of higher, the equispaced libMesh nodes no longer are a subset of the GLL ndoes.

Constructor & Destructor Documentation

◆ NekRSMesh() [1/2]

NekRSMesh::NekRSMesh ( const InputParameters &  parameters)

◆ NekRSMesh() [2/2]

NekRSMesh::NekRSMesh ( const NekRSMesh )
default

Member Function Documentation

◆ addElems()

virtual void NekRSMesh::addElems ( )
virtual

Add all the elements in the mesh to the MOOSE data structures.

◆ boundary()

const std::vector<int>* NekRSMesh::boundary ( ) const
inline

Get the boundary ID for which nekRS and MOOSE are coupled

Returns
boundary ID

◆ boundary_id()

int NekRSMesh::boundary_id ( const int  elem_id,
const int  face_id 
)
protected

Sideset ID corresponding to a given volume element with give local face ID

Parameters
[in]elem_idelement local rank ID
[in]face_idelement-local face ID
Returns
sideset ID (-1 means not one a boundary)

◆ boundaryCoupling()

const NekBoundaryCoupling& NekRSMesh::boundaryCoupling ( ) const
inline

Get the boundary coupling data structure

Returns
boundary coupling data structure

◆ boundaryElem()

Elem* NekRSMesh::boundaryElem ( ) const

Create a new element for a boundary mesh.

◆ boundaryElemProcessorID()

int NekRSMesh::boundaryElemProcessorID ( const int  elem_id)

Processor id (rank) owning the given boundary element

Returns
processor id

◆ boundaryNodeIndex()

int NekRSMesh::boundaryNodeIndex ( const int  gll_index) const
inline

Get the libMesh node index from nekRS's GLL index ordering

Parameters
[in]gll_indexnekRS GLL index
Returns
node index

◆ buildDummyMesh()

virtual void NekRSMesh::buildDummyMesh ( )
virtual

If running NekRS in JIT mode, we still need to make a mesh based on requirements in MOOSE, so we just make a dummy mesh of a single Quad4 element

◆ buildMesh()

virtual void NekRSMesh::buildMesh ( )
overridevirtual

◆ cornerIndices()

std::vector<std::vector<int> > NekRSMesh::cornerIndices ( ) const
inline

Get the corner indices for the GLL points to be used in the mesh mirror

Returns
mapping of mesh mirror nodes to GLL points

◆ exactMirror()

bool NekRSMesh::exactMirror ( ) const
inline

Whether the mesh mirror is an exact representation of the NekRS mesh

Parameters
returnwhether mesh mirror is exact

◆ extractSurfaceMesh()

virtual void NekRSMesh::extractSurfaceMesh ( )
virtual

For the case of surface coupling only (i.e. no volume coupling), we create a surface mesh for the elements on the specified boundary IDs

◆ extractVolumeMesh()

virtual void NekRSMesh::extractVolumeMesh ( )
virtual

For the case of volume coupling only (i.e. no surface coupling), we create a volume mesh for all volume elements

◆ facesOnBoundary()

int NekRSMesh::facesOnBoundary ( const int  elem_id) const

Get the number of faces of this global element that are on a coupling boundary

Parameters
[in]elem_idglobal element ID
Returns
number of faces on a coupling boundary

◆ faceVertices()

void NekRSMesh::faceVertices ( )
protected

Get the vertices defining the surface mesh interpolation from the stored coupling information and store in _x, _y, and _z

◆ initializeMeshParams()

void NekRSMesh::initializeMeshParams ( )
protected

Initialize members for the mesh and determine the GLL-to-node mapping.

◆ initializePreviousDisplacements()

void NekRSMesh::initializePreviousDisplacements ( )

Initialize previous displacement values to zero for boundary mirror-based moving mesh problems.

◆ nBuildPerSurfaceElem()

int NekRSMesh::nBuildPerSurfaceElem ( ) const
inline

Number of elements to build for each NekRS surface element

Returns
number of MOOSE elements per NekRS surface element

◆ nBuildPerVolumeElem()

int NekRSMesh::nBuildPerVolumeElem ( ) const
inline

Number of elements to build for each NekRS volume element

Returns
number of MOOSE elements per NekRS volume element

◆ nek_initial_x()

const std::vector<double>& NekRSMesh::nek_initial_x ( ) const
inline

Get the initial mesh x coordinates

Returns
initial mesh x coordinates

◆ nek_initial_y()

const std::vector<double>& NekRSMesh::nek_initial_y ( ) const
inline

Get the initial mesh y coordinates

Returns
initial mesh y coordinates

◆ nek_initial_z()

const std::vector<double>& NekRSMesh::nek_initial_z ( ) const
inline

Get the initial mesh z coordinates

Returns
initial mesh z coordinates

◆ nekNumQuadraturePoints1D()

int NekRSMesh::nekNumQuadraturePoints1D ( ) const

Get the number of quadrature points per coordiante direction in nekRS's mesh

Returns
number of quadrature points per coordinate direction

◆ nekNumSurfaceElems()

const int& NekRSMesh::nekNumSurfaceElems ( ) const
inline

Get the total number of surface elements in nekRS's mesh

Returns
number of surface elements

◆ nekNumVolumeElems()

const int& NekRSMesh::nekNumVolumeElems ( ) const
inline

Get the total number of volume elements in nekRS's mesh

Returns
number of volume elements

◆ nekPolynomialOrder()

int NekRSMesh::nekPolynomialOrder ( ) const
inline

NekRS mesh polynomial order

Returns
NekRS polynomial order

◆ nMoosePerNek()

int NekRSMesh::nMoosePerNek ( ) const
inline

Get the number of MOOSE elements we build for each NekRS element

Returns
MOOSE elements built for each NekRS element

◆ nodeIndex()

int NekRSMesh::nodeIndex ( const int  gll_index) const
inline

Get the libMesh node index from nekRS's GLL index ordering.

This function is used to perform the data transfer routines in NekRSProblem agnostic of whether we have surface or volume coupling.

Parameters
[in]gll_indexnekRS GLL index
Returns
node index

◆ numElems()

const int& NekRSMesh::numElems ( ) const
inline

Get the number of NekRS elements we rebuild in the MOOSE mesh.

This function is used to perform the data transfer routines in NekRSProblem agnostic of whether we have surface or volume coupling. return number of elements

◆ numQuadraturePoints1D()

int NekRSMesh::numQuadraturePoints1D ( ) const

Get the number of quadrature points per coordinate direction in MOOSE's representation of nekRS's mesh

Returns
number of quadrature points per coordinate direction

◆ numSurfaceElems()

const int& NekRSMesh::numSurfaceElems ( ) const
inline

Get the number of surface elements in MOOSE's representation of nekRS's mesh

Returns
number of surface elements

◆ numVerticesPerElem()

const int& NekRSMesh::numVerticesPerElem ( ) const
inline

Get the number of vertices per element in MOOSE's representation of nekRS's mesh.

This function is used to perform the data transfer routines in NekRSProblem agnostic of whether we have surface or volume coupling. return number of vertices per element

◆ numVerticesPerSurface()

const int& NekRSMesh::numVerticesPerSurface ( ) const
inline

Get the number of vertices per surface element in MOOSE's representation of nekRS's mesh

Returns
number of vertices per surface

◆ numVerticesPerVolume()

const int& NekRSMesh::numVerticesPerVolume ( ) const
inline

Get the number of vertices per volume element in MOOSE's representation of nekRS's mesh

Returns
number of vertices per volume

◆ numVolumeElems()

const int& NekRSMesh::numVolumeElems ( ) const
inline

Get the number of volume elements in MOOSE's representation of nekRS's mesh

Returns
number of volume elements

◆ operator=()

NekRSMesh& NekRSMesh::operator= ( const NekRSMesh other_mesh)
delete

◆ order()

const order::NekOrderEnum& NekRSMesh::order ( ) const
inline

Get the order of the surface mesh; note that this is zero-indexed, so 0 = first, 1 = second

Returns
order

◆ prev_disp_x()

std::vector<double>& NekRSMesh::prev_disp_x ( )
inline

Get the previous x displacement

Returns
previous x displacement values

◆ prev_disp_y()

std::vector<double>& NekRSMesh::prev_disp_y ( )
inline

Get the previous y displacement

Returns
previous y displacement values

◆ prev_disp_z()

std::vector<double>& NekRSMesh::prev_disp_z ( )
inline

Get the previous z displacement

Returns
previous z displacement values

◆ printMeshInfo()

virtual void NekRSMesh::printMeshInfo ( ) const
virtual

Print diagnostic information related to the mesh.

◆ safeClone()

virtual std::unique_ptr<MooseMesh> NekRSMesh::safeClone ( ) const
overridevirtual

◆ saveInitialVolMesh()

void NekRSMesh::saveInitialVolMesh ( )

Save the initial volumetric mesh for volume mirror-based moving mesh problems.

◆ scaling()

const Real& NekRSMesh::scaling ( ) const
inline

Get the scaling factor applied to the nekRS mesh

Returns
scaling factor

◆ storeBoundaryCoupling()

void NekRSMesh::storeBoundaryCoupling ( )
protected

Store the rank-local element and rank ownership for boundary coupling; this loops over the NekRS mesh and fetches relevant information on the boundaries

◆ storeVolumeCoupling()

void NekRSMesh::storeVolumeCoupling ( )
protected

Store the rank-local element and rank ownership for volume coupling.

◆ updateDisplacement()

void NekRSMesh::updateDisplacement ( const int  e,
const double *  src,
const field::NekWriteEnum  field 
)

Copy a new boundary displacement value for a given element's face

Parameters
[in]srcthe source displacement at element e and face f
[in]ethe element to which src belongs
[in]fieldthe displacement field we are updating

◆ validParams()

static InputParameters NekRSMesh::validParams ( )
static

◆ volume()

const bool& NekRSMesh::volume ( ) const
inline

Get whether the mesh permits volume-based coupling

Returns
whether mesh contains volume elements

◆ volumeCoupling()

const NekVolumeCoupling& NekRSMesh::volumeCoupling ( ) const
inline

Get the volume coupling data structure

Returns
volume coupling data structure

◆ volumeElem()

Elem* NekRSMesh::volumeElem ( ) const

Create a new element for a volume mesh.

◆ volumeElemProcessorID()

int NekRSMesh::volumeElemProcessorID ( const int  elem_id)

Processor id (rank) owning the given volume element

Returns
processor id

◆ volumeNodeIndex()

int NekRSMesh::volumeNodeIndex ( const int  gll_index) const
inline

Get the libMesh node index from nekRS's GLL index ordering

Parameters
[in]gll_indexnekRS GLL index
Returns
node index

◆ volumeVertices()

void NekRSMesh::volumeVertices ( )
protected

Get the vertices defining the volume mesh interpolation from the stored coupling information and store in _x, _y, and _z

Member Data Documentation

◆ _bnd_node_index

std::vector<int> NekRSMesh::_bnd_node_index
protected

Mapping of boundary GLL indices to MooseMesh node indices.

In nekRS, the GLL points are ordered by \(x\), \(y\), and \(z\) coordinates, but in order to construct sensible elements in Moose, we need to reorder these points so that they match a libMesh-friendly node ordering. Without such a mapping, we would construct triangles with zero/negative Jacobians instead of quad elements. By indexing in the GLL index, this returns the node index.

◆ _boundary

const std::vector<int>* NekRSMesh::_boundary
protected

Boundary ID(s) through which to couple Nek to MOOSE.

◆ _boundary_coupling

NekBoundaryCoupling NekRSMesh::_boundary_coupling
protected

Data structure holding mapping information for boundary coupling.

◆ _corner_indices

std::vector<std::vector<int> > NekRSMesh::_corner_indices
protected

Corner indices for GLL points of mesh mirror elements.

◆ _elem_processor_id

int(NekRSMesh::* NekRSMesh::_elem_processor_id) (const int elem_id)
protected

Function returning the processor id which should own each element.

◆ _exact

const bool& NekRSMesh::_exact
protected

Whether the NekRS mesh mirror is an exact replica of the NekRS mesh. If false (the default), then we build one MOOSE element for each NekRS element, and the order of the MOOSE element is selected with 'order'. If true, then we instead build one MOOSE element for each "first-order element" within each NekRS high-order spectral element. In other words, if the NekRS mesh is polynomial order 7, then we would build 7^2 MOOSE surface elements for each NekRS surface element, and 7^3 MOOSE volume elements for each NekRS volume element. The order of these elements will be first order.

◆ _fluid_block_id

const unsigned int& NekRSMesh::_fluid_block_id
protected

Block ID for the fluid portion of the mesh mirror.

◆ _initial_x

std::vector<double> NekRSMesh::_initial_x
protected

\(x\) coordinates of the initial GLL points, for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.

◆ _initial_y

std::vector<double> NekRSMesh::_initial_y
protected

\(y\) coordinates of the initial GLL points, for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.

◆ _initial_z

std::vector<double> NekRSMesh::_initial_z
protected

\(z\) coordinates of the initial GLL points, for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.

◆ _n_build_per_surface_elem

int NekRSMesh::_n_build_per_surface_elem
protected

Number of MOOSE surface elements to build per NekRS surface element.

◆ _n_build_per_volume_elem

int NekRSMesh::_n_build_per_volume_elem
protected

Number of MOOSE volume elements to build per NekRS volume element.

◆ _n_elems

int NekRSMesh::_n_elems
protected

Number of NekRS elements we build in the MooseMesh. The total number of elements in the mesh mirror is _n_elems * _n_moose_per_nek.

◆ _n_moose_per_nek

int NekRSMesh::_n_moose_per_nek
protected

Number of MOOSE elements corresponding to each NekRS element, which depends on whether building a boundary/volume mesh

◆ _n_surface_elems

int NekRSMesh::_n_surface_elems
protected

Number of NekRS surface elements in MooseMesh. The total number of surface elements in the mesh mirror is _n_surface_elems * _n_build_per_surface_elem.

◆ _n_vertices_per_elem

int NekRSMesh::_n_vertices_per_elem
protected

Number of vertices per element, which depends on whether building a boundary/volume mesh.

◆ _n_vertices_per_surface

int NekRSMesh::_n_vertices_per_surface
protected

Number of vertices per surface.

◆ _n_vertices_per_volume

int NekRSMesh::_n_vertices_per_volume
protected

Number of vertices per volume element.

◆ _n_volume_elems

int NekRSMesh::_n_volume_elems
protected

Number of NekRS volume elements in MooseMesh. The total number of volume elements in the mesh mirror is _n_volume_elems * _n_build_per_volume_elem.

◆ _nek_internal_mesh

mesh_t* NekRSMesh::_nek_internal_mesh = nullptr
protected

Pointer to NekRS's internal mesh data structure.

◆ _nek_n_surface_elems

int NekRSMesh::_nek_n_surface_elems
protected

Total number of surface elements in the nekRS problem.

◆ _nek_n_volume_elems

int NekRSMesh::_nek_n_volume_elems
protected

Total number of volume elements in the nekRS problem.

◆ _nek_polynomial_order

int NekRSMesh::_nek_polynomial_order
protected

Order of the nekRS solution.

◆ _new_elem

Elem*(NekRSMesh::* NekRSMesh::_new_elem) () const
protected

Function pointer to the type of new element to add.

◆ _node_index

std::vector<int>* NekRSMesh::_node_index
protected

Mapping of GLL nodes to libMesh nodes, which depends on whether building a boundary/volume mesh.

◆ _order

const order::NekOrderEnum NekRSMesh::_order
protected

Order of the surface interpolation between nekRS and MOOSE.

Order of the interpolation to be performed between nekRS and MOOSE; options = FIRST, SECOND. For a first-order interpolation, nekRS's solution is interpolated onto a first-order surface mesh (i.e. Quad4), while for a second-order interpolation, nekRS's solution is interpolated onto a second-order surface mesh (i.e. Quad9). Note that this is zero-indexed so that an integer value of 0 = first-order, 1 = second-order, etc.

◆ _phase

std::vector<int> NekRSMesh::_phase
protected

"Phase" for each element (fluid = 0, solid = 1)

TODO: This could be improved so that we don't need to collect this information, but would require rewriting of a lot of the code used to assemble the members in NekVolumeCoupling.

◆ _prev_disp_x

std::vector<double> NekRSMesh::_prev_disp_x
protected

\(x\), \(y\), \(z\) displacements of the boundary mesh mirror for calculating displacement, in this rank

◆ _prev_disp_y

std::vector<double> NekRSMesh::_prev_disp_y
protected

\(x\), \(y\), \(z\) displacements of the boundary mesh mirror for calculating displacement, in this rank

◆ _prev_disp_z

std::vector<double> NekRSMesh::_prev_disp_z
protected

\(x\), \(y\), \(z\) displacements of the boundary mesh mirror for calculating displacement, in this rank

◆ _scaling

const Real& NekRSMesh::_scaling
protected

Spatial scaling factor to apply to the mesh.

nekRS is dimension agnostic - depending on the values used for the material properties, the units of the mesh are arbitrary. Other apps that nekRS might be coupled to could be in different units - to allow each app to use the units that it wants, we can simply scale the NekRSMesh by a constant factor. This will also adjust the heat flux coming in to nekRS by an appropriate factor. For instance, if nekRS solves a problem in units of meters, but a BISON solution is done on a mesh in units of centimeters, this scaling factor should be set to

  1. Note that other postprocessors will still be calculated on the nekRS mesh, which will be in whatever units nekRS is internally using.

◆ _side_index

std::vector<int> NekRSMesh::_side_index
protected

Mapping of side indices to libMesh side indices.

nekRS uses its own side mapping that differs from that assumed in libMesh. In order to assign the correct sideset IDs to the MooseMesh, we need to know the mapping between these different conventions. By indexing in the nekRS side index, this returns the libMesh side index.

◆ _solid_block_id

const unsigned int& NekRSMesh::_solid_block_id
protected

Block ID for the solid portion of the mesh mirror.

◆ _vol_node_index

std::vector<int> NekRSMesh::_vol_node_index
protected

Mapping of volume GLL indices to MooseMesh node indices.

In nekRS, the GLL points are ordered by \(x\), \(y\), and \(z\) coordinates, but in order to construct sensible elements in Moose, we need to reorder these points so that they match a libMesh-friendly node ordering. Without such a mapping, we would construct triangles with zero/negative Jacobians instead of hex elements. By indexing in the GLL index, this returns the node index.

◆ _volume

const bool& NekRSMesh::_volume
protected

Whether nekRS is coupled through volumes to MOOSE.

Unlike the case with _boundary, nekRS has no concept of volume/block IDs, so we cannot have the user provide a vector of volumes that they want to construct, so the best we can do is use a boolean here to turn on/off the volume-based coupling for the entire mesh.

◆ _volume_coupling

NekVolumeCoupling NekRSMesh::_volume_coupling
protected

Data structure holding mapping information for volume coupling.

◆ _x

std::vector<double> NekRSMesh::_x
protected

\(x\) coordinates of the current GLL points (which can move in time), for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.

◆ _y

std::vector<double> NekRSMesh::_y
protected

\(y\) coordinates of the current GLL points (which can move in time), for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.

◆ _z

std::vector<double> NekRSMesh::_z
protected

\(z\) coordinates of the current GLL points (which can move in time), for this rank

This is ordered according to nekRS's internal geometry layout, and is indexed first by the element and then by the node.


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