Cardinal
Classes | Functions | Variables
nekrs Namespace Reference

Cardinal-specific nekRS API. More...

Classes

struct  characteristicScales
 Characteristic scales assumed in nekRS if using a non-dimensional solution. More...
 
struct  usrwrkIndices
 

Functions

std::vector< std::vector< int > > cornerGLLIndices (const int &n, const bool exact)
 Get the corner indices for the GLL points in a NekRS element on a single face. More...
 
std::vector< std::vector< int > > nestedElementsOnFace (const int &n)
 
void initializeHostMeshParameters ()
 Allocate memory for the host mesh parameters. More...
 
void updateHostMeshParameters ()
 Update the mesh parameters on host. More...
 
dfloat * getSgeo ()
 
dfloat * getVgeo ()
 
void setAbsoluteTol (double tol)
 
void setRelativeTol (double tol)
 
void setNekSetupTime (const double &time)
 
double getNekSetupTime ()
 
void setStartTime (const double &start)
 
bool isInitialized ()
 
void write_usrwrk_field_file (const int &slot, const std::string &prefix, const dfloat &time, const int &step, const bool &write_coords)
 
void write_field_file (const std::string &prefix, const dfloat time, const int &step)
 
void buildOnly (int buildOnly)
 
int buildOnly ()
 
void interpolateVolumeHex3D (const double *I, double *x, int N, double *Ix, int M)
 
bool hasCHT ()
 
bool hasMovingMesh ()
 
bool hasVariableDt ()
 
bool hasBlendingSolver ()
 
bool hasUserMeshSolver ()
 
bool endControlElapsedTime ()
 
bool endControlTime ()
 
bool endControlNumSteps ()
 
int scalarFieldOffset ()
 
int velocityFieldOffset ()
 
mesh_t * entireMesh ()
 
mesh_t * flowMesh ()
 
mesh_t * temperatureMesh ()
 
mesh_t * getMesh (const nek_mesh::NekMeshEnum pp_mesh)
 
int commRank ()
 
int commSize ()
 
bool hasTemperatureVariable ()
 
bool hasTemperatureSolve ()
 
bool hasScalarVariable (int scalarId)
 
bool hasHeatSourceKernel ()
 
bool scratchAvailable ()
 
void initializeScratch (const unsigned int &n_slots)
 
void freeScratch ()
 Free the scratch space. More...
 
double viscosity ()
 
double Pr ()
 
void copyDeformationToDevice ()
 Copy the deformation from host to device. More...
 
template<typename T >
void allgatherv (const std::vector< int > &base_counts, const T *input, T *output, const int multiplier=1)
 
void displacementAndCounts (const std::vector< int > &base_counts, int *counts, int *displacement, const int multiplier)
 
void interpolationMatrix (double *I, int starting_points, int ending_points)
 
void interpolateSurfaceFaceHex3D (double *scratch, const double *I, double *x, int N, double *Ix, int M)
 
Point centroidFace (int local_elem_id, int local_face_id)
 
Point centroid (int local_elem_id)
 
Point gllPoint (int local_elem_id, int local_node_id)
 
Point gllPointFace (int local_elem_id, int local_face_id, int local_node_id)
 
std::vector< double > usrwrkSideIntegral (const unsigned int &slot, const std::vector< int > &boundary, const nek_mesh::NekMeshEnum pp_mesh)
 
double usrwrkVolumeIntegral (const unsigned int &slot, const nek_mesh::NekMeshEnum pp_mesh)
 
void scaleUsrwrk (const unsigned int &slot, const dfloat &value)
 
bool normalizeFluxBySideset (const NekBoundaryCoupling &nek_boundary_coupling, const std::vector< int > &boundary, const std::vector< double > &moose_integral, std::vector< double > &nek_integral, double &normalized_nek_integral)
 
bool normalizeFlux (const NekBoundaryCoupling &nek_boundary_coupling, const std::vector< int > &boundary, const double moose_integral, double nek_integral, double &normalized_nek_integral)
 
double area (const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
 
double sideIntegral (const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
 
double volume (const nek_mesh::NekMeshEnum pp_mesh)
 
void dimensionalizeVolume (double &integral)
 
void dimensionalizeArea (double &integral)
 
void dimensionalizeVolumeIntegral (const field::NekFieldEnum &integrand, const Real &volume, double &integral)
 
void dimensionalizeSideIntegral (const field::NekFieldEnum &integrand, const Real &area, double &integral)
 
void dimensionalizeSideIntegral (const field::NekFieldEnum &integrand, const std::vector< int > &boundary_id, double &integral, const nek_mesh::NekMeshEnum pp_mesh)
 
double volumeIntegral (const field::NekFieldEnum &integrand, const double &volume, const nek_mesh::NekMeshEnum pp_mesh)
 
double massFlowrate (const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
 
double sideMassFluxWeightedIntegral (const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
 
double pressureSurfaceForce (const std::vector< int > &boundary_id, const Point &direction, const nek_mesh::NekMeshEnum pp_mesh)
 
double heatFluxIntegral (const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
 
void limitTemperature (const double *min_T, const double *max_T)
 
void gradient (const int offset, const double *f, double *grad_f, const nek_mesh::NekMeshEnum pp_mesh)
 
double volumeExtremeValue (const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
 
double sideExtremeValue (const std::vector< int > &boundary_id, const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
 
int Nfaces ()
 
bool isHeatFluxBoundary (const int boundary)
 
bool isMovingMeshBoundary (const int boundary)
 
bool isTemperatureBoundary (const int boundary)
 
const std::string temperatureBoundaryType (const int boundary)
 
int polynomialOrder ()
 
int Nelements ()
 
int dim ()
 
int NfaceVertices ()
 Number of vertices required to define an element face Vertices refer to the points required to place the "corners" of an element face, and not the quadrature points. For instance, for hexahedral elements, the number of vertices per face is 4 regardless of the polynomial order. More...
 
int NboundaryFaces ()
 
int NboundaryID ()
 
bool validBoundaryIDs (const std::vector< int > &boundary_id, int &first_invalid_id, int &n_boundaries)
 
void storeBoundaryCoupling (const std::vector< int > &boundary_id, int &N)
 
double scalar01 (const int id)
 Get the scalar01 solution at given GLL index. More...
 
double scalar02 (const int id)
 Get the scalar02 solution at given GLL index. More...
 
double scalar03 (const int id)
 Get the scalar03 solution at given GLL index. More...
 
double temperature (const int id)
 Get the temperature solution at given GLL index. More...
 
double pressure (const int id)
 
double unity (const int id)
 
double velocity_x (const int id)
 
double velocity_y (const int id)
 
double velocity_z (const int id)
 
double velocity (const int id)
 
void flux (const int id, const dfloat value)
 
void heat_source (const int id, const dfloat value)
 
void x_displacement (const int id, const dfloat value)
 
void y_displacement (const int id, const dfloat value)
 
void z_displacement (const int id, const dfloat value)
 
void initializeDimensionalScales (const double U_ref, const double T_ref, const double dT_ref, const double L_ref, const double rho_ref, const double Cp_ref)
 
void dimensionalize (const field::NekFieldEnum &field, double &value)
 Dimensionalize a field by multiplying the nondimensional form by the reference. More...
 
double referenceFlux ()
 
double referenceSource ()
 
double referenceLength ()
 
double referenceArea ()
 
double referenceVolume ()
 
template<typename T >
MPI_Datatype resolveType ()
 

Variables

double(*)(int) solutionPointer (const field::NekFieldEnum &field)
 

Detailed Description

Cardinal-specific nekRS API.

nekRS ships with a rudimentary API in their nekrs namespace, but we need additional functionality from within Cardinal. Many of these functions are quite basic and could eventually be ported back into nekRS itself.

Function Documentation

◆ allgatherv()

template<typename T >
void nekrs::allgatherv ( const std::vector< int > &  base_counts,
const T *  input,
T *  output,
const int  multiplier 
)

Helper function for MPI_Allgatherv of results in NekRS

Parameters
[in]base_countsonce multiplied by 'multiplier', the number of counts on each rank
[in]inputrank-local data
[out]outputcollected result
[in]multiplierconstant multiplier to set on each count indicator

◆ area()

double nekrs::area ( const std::vector< int > &  boundary_id,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the area of a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to perform the integral
[in]pp_meshwhich NekRS mesh to operate on
Returns
area integral

◆ buildOnly() [1/2]

int nekrs::buildOnly ( )

Whether NekRS was run in JIT build-only mode

Returns
whether NekRS was run in build-only mode

◆ buildOnly() [2/2]

void nekrs::buildOnly ( int  buildOnly)

Indicate whether NekRS was run in build-only mode (this doesn't actually cause NekRS to run in build-only mode, but only provides an interface to this information elsewhere).

Parameters
[in]buildOnlywhether NekRS is to be run in build-only mode

◆ centroid()

Point nekrs::centroid ( int  local_elem_id)

Compute the centroid given a local element ID (NOTE: returns in dimensional form)

Parameters
[in]local_elem_idlocal element ID on this rank
Returns
centroid

◆ centroidFace()

Point nekrs::centroidFace ( int  local_elem_id,
int  local_face_id 
)

Compute the face centroid given a local element ID and face ID (NOTE: returns in dimensional form)

Parameters
[in]local_elem_idlocal element ID on this rank
[in]local_face_idlocal face ID on the element
Returns
centroid

◆ commRank()

int nekrs::commRank ( )

Get the process rank

Returns
process rank

◆ commSize()

int nekrs::commSize ( )

Get the communicator size

Returns
communicator size

◆ copyDeformationToDevice()

void nekrs::copyDeformationToDevice ( )

Copy the deformation from host to device.

◆ cornerGLLIndices()

std::vector< std::vector< int > > nekrs::cornerGLLIndices ( const int &  n,
const bool  exact 
)

Get the corner indices for the GLL points in a NekRS element on a single face.

Support the NekRS mesh face has nodes

6 – 7 – 8 | | | 3 – 4 – 5 | | | 0 – 1 – 2

If 'exact' is false, this method returns {0, 2, 6, 8}. If 'exact' is true, this method returns:

  • element 0: {0, 1, 3, 4}
  • element 1: {1, 2, 4, 5}
  • element 2: {3, 4, 6, 7}
  • element 3: {4, 5, 7, 8}
Parameters
[in]norder of mesh (1 = first-order, 2 = second-order)
[in]exactwhether the MOOSE elements will exactly represent the NekRS mesh
Returns
GLL indices, indexed first by element to build, then second by nodes on that element

◆ dim()

int nekrs::dim ( )

Mesh dimension

Returns
mesh dimension

◆ dimensionalize()

void nekrs::dimensionalize ( const field::NekFieldEnum field,
double &  value 
)

Dimensionalize a field by multiplying the nondimensional form by the reference.

This routine dimensionalizes a nondimensional term by multiplying the non-dimensional form by a scalar, i.e. \(f^\dagger=\frac{f}{f_ref}\), where \(f^\dagger\) is the nondimensional form and \(f_{ref}\) is a reference scale with form particular to the interpretation of the field. Note that for temperature in particular, there are still additional steps to dimensionalize, because we do not define a nondimensional temperature simply as \(T^\dagger=\frac{T}{\Delta T_{ref}}\). But, this function just treats the characteristic scale that would appear in the denominator.

Parameters
[in]fieldphysical interpretation of value to dimensionalize
[out]valuevalue to dimensionalize

◆ dimensionalizeArea()

void nekrs::dimensionalizeArea ( double &  integral)

Dimensionalize an area

Parameters
[in]integralintegral to dimensionalize

◆ dimensionalizeSideIntegral() [1/2]

void nekrs::dimensionalizeSideIntegral ( const field::NekFieldEnum integrand,
const Real &  area,
double &  integral 
)

Dimensionalize a given integral of f over a side, i.e. fdS

Parameters
[in]integrandfield to dimensionalize
[in]areaarea of the boundary
[in]integralintegral to dimensionalize

◆ dimensionalizeSideIntegral() [2/2]

void nekrs::dimensionalizeSideIntegral ( const field::NekFieldEnum integrand,
const std::vector< int > &  boundary_id,
double &  integral,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Dimensionalize a given integral of f over a side, i.e. fdS

Parameters
[in]integrandfield to dimensionalize
[in]boundary_idboundary IDs for the integral
[in]integralintegral to dimensionalize
[in]pp_meshwhich NekRS mesh to operate on

◆ dimensionalizeVolume()

void nekrs::dimensionalizeVolume ( double &  integral)

Dimensionalize a volume

Parameters
[in]integralintegral to dimensionalize

◆ dimensionalizeVolumeIntegral()

void nekrs::dimensionalizeVolumeIntegral ( const field::NekFieldEnum integrand,
const Real &  volume,
double &  integral 
)

Dimensionalize a given integral of f over volume, i.e. fdV

Parameters
[in]integrandfield to dimensionalize
[in]volumevolume of the domain (only used for dimensionalizing temperature)
[in]integralintegral to dimensionalize

◆ displacementAndCounts()

void nekrs::displacementAndCounts ( const std::vector< int > &  base_counts,
int *  counts,
int *  displacement,
const int  multiplier 
)

Determine the receiving counts and displacements for all gather routines

Parameters
[in]base_countsunit-wise receiving counts for each process
[out]countsreceiving counts from each process
[out]displacementdisplacement for each process's counts
[in]multiplieroptional multiplier on the face-based data

◆ endControlElapsedTime()

bool nekrs::endControlElapsedTime ( )

Whether nekRS's input file intends to terminate the simulation based on a wall time

Returns
whether a wall time is used in nekRS to end the simulation

◆ endControlNumSteps()

bool nekrs::endControlNumSteps ( )

Whether nekRS's input file intends to terminate the simulation based on a number of steps

Returns
whether a time step interval is used in nekRS to end the simulation

◆ endControlTime()

bool nekrs::endControlTime ( )

Whether nekRS's input file intends to terminate the simulation based on an end time

Returns
whether an end time is used in nekRS to end the simulation

◆ entireMesh()

mesh_t* nekrs::entireMesh ( )

Get the "entire" NekRS mesh. For cases with a temperature scalar, this returns nrs->meshT, which will cover both the fluid and solid regions if CHT is present. For flow-only cases, this will return the flow mesh.

Returns
entire NekRS mesh

◆ flowMesh()

mesh_t* nekrs::flowMesh ( )

Get the mesh for the flow solve

Returns
flow mesh

◆ flux()

void nekrs::flux ( const int  id,
const dfloat  value 
)

Write a value into the user scratch space that holds the flux

Parameters
[in]idindex
[in]valuevalue to write

◆ freeScratch()

void nekrs::freeScratch ( )

Free the scratch space.

◆ getMesh()

mesh_t* nekrs::getMesh ( const nek_mesh::NekMeshEnum  pp_mesh)

Get the mesh to act on

Parameters
[in]pp_meshwhich NekRS mesh to operate on
Returns
mesh to act on

◆ getNekSetupTime()

double nekrs::getNekSetupTime ( )

Get time spent on initialization

Returns
time spent on initialization

◆ getSgeo()

dfloat* nekrs::getSgeo ( )

◆ getVgeo()

dfloat* nekrs::getVgeo ( )

◆ gllPoint()

Point nekrs::gllPoint ( int  local_elem_id,
int  local_node_id 
)

Get the coordinate given a local element ID and local node ID (NOTE: returns in dimensional form)

Parameters
[in]local_elem_idlocal element ID on this rank
[in]local_node_idlocal node ID on this element
Returns
point

◆ gllPointFace()

Point nekrs::gllPointFace ( int  local_elem_id,
int  local_face_id,
int  local_node_id 
)

Get the coordinate given a local element ID, a local face ID, and local node ID (NOTE: returns in dimensional form)

Parameters
[in]local_elem_idlocal element ID on this rank
[in]local_face_idlocal face ID on this element
[in]local_node_idlocal node ID on this element
Returns
point

◆ gradient()

void nekrs::gradient ( const int  offset,
const double *  f,
double *  grad_f,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the gradient of a volume field

Parameters
[in]offsetin the gradient field for each component (grad_x, grad_y, or grad_z)
[in]ffield to compute the gradient of
[in]pp_meshwhich NekRS mesh to operate on
[out]grad_fgradient of field

◆ hasBlendingSolver()

bool nekrs::hasBlendingSolver ( )

Whether nekRS's input file has the blending mesh solver

Returns
whether nekRS's input file has a non-user [MESH] solver

◆ hasCHT()

bool nekrs::hasCHT ( )

Whether nekRS's input file has CHT

Returns
whether nekRS input files model CHT

◆ hasHeatSourceKernel()

bool nekrs::hasHeatSourceKernel ( )

Whether nekRS contains an OCCA kernel to apply a source to the passive scalar equations

Returns
whether nekRS has an OCCA kernel for apply a passive scalar source

◆ hasMovingMesh()

bool nekrs::hasMovingMesh ( )

Whether nekRS's input file indicates a moving mesh

Returns
whether nekRS's input file indicates a moving mesh

◆ hasScalarVariable()

bool nekrs::hasScalarVariable ( int  scalarId)

Whether nekRS's input file indicates that the problem has a scalar0(scalarId) variable

Parameters
[in]scalarIdscalar number, i.e. for scalar03 scalarId=3
Returns
whether the nekRS problem includes the scalar0(scalarId) variable

◆ hasTemperatureSolve()

bool nekrs::hasTemperatureSolve ( )

Whether nekRS actually solves for temperature (as opposed to setting its solver to 'none')

Returns
whether nekRS will solve for temperature

◆ hasTemperatureVariable()

bool nekrs::hasTemperatureVariable ( )

Whether nekRS's input file indicates that the problem has a temperature variable

Returns
whether the nekRS problem includes a temperature variable

◆ hasUserMeshSolver()

bool nekrs::hasUserMeshSolver ( )

Whether nekRS's input file has the user mesh solver

Returns
whether nekRS's input file has [MESH] solver = user

◆ hasVariableDt()

bool nekrs::hasVariableDt ( )

Whether nekRS's input file indicates a variable time stepping scheme

Returns
whether nekRS's input file indicates a variable time stepping

◆ heat_source()

void nekrs::heat_source ( const int  id,
const dfloat  value 
)

Write a value into the user scratch space that holds the volumetric heat source

Parameters
[in]idindex
[in]valuevalue to write

◆ heatFluxIntegral()

double nekrs::heatFluxIntegral ( const std::vector< int > &  boundary_id,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the heat flux over a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to perform the integral
[in]pp_meshwhich NekRS mesh to operate on
Returns
heat flux area integral

◆ initializeDimensionalScales()

void nekrs::initializeDimensionalScales ( const double  U_ref,
const double  T_ref,
const double  dT_ref,
const double  L_ref,
const double  rho_ref,
const double  Cp_ref 
)

Initialize the characteristic scales for a nondimesional solution

Parameters
[in]U_refreference velocity
[in]T_refreference temperature
[in]dT_refreference temperature range
[in]L_refreference length scale
[in]rho_refreference density
[in]Cp_refreference heat capacity

◆ initializeHostMeshParameters()

void nekrs::initializeHostMeshParameters ( )

Allocate memory for the host mesh parameters.

◆ initializeScratch()

void nekrs::initializeScratch ( const unsigned int &  n_slots)

Initialize scratch space for data to get sent into NekRS

Parameters
[in]n_slotsnumber of slots (for volume arrays) to allocate

◆ interpolateSurfaceFaceHex3D()

void nekrs::interpolateSurfaceFaceHex3D ( double *  scratch,
const double *  I,
double *  x,
int  N,
double *  Ix,
int  M 
)

Interpolate face data onto a new set of points

Parameters
[in]scratchavailable scratch space for the calculation
[in]Iinterpolation matrix
[in]xface data to be interpolated
[in]Nnumber of points in 1-D to be interpolated
[out]Ixinterpolated data
[in]Mresulting number of interpolated points in 1-D

◆ interpolateVolumeHex3D()

void nekrs::interpolateVolumeHex3D ( const double *  I,
double *  x,
int  N,
double *  Ix,
int  M 
)

Interpolate a volume between NekRS's GLL points and a given-order receiving/sending mesh

◆ interpolationMatrix()

void nekrs::interpolationMatrix ( double *  I,
int  starting_points,
int  ending_points 
)

Form the 2-D interpolation matrix from a starting GLL quadrature rule to an ending GLL quadrature rule.

Parameters
[out]Iinterpolation matrix
[in]starting_pointsnumber of points in the source quadrature rule
[in]ending_pointsnumber of points in the end quadrature rule

◆ isHeatFluxBoundary()

bool nekrs::isHeatFluxBoundary ( const int  boundary)

Whether the specific boundary is a flux boundary

Parameters
[in]boundaryboundary ID
Returns
whether boundary is a flux boundary

◆ isInitialized()

bool nekrs::isInitialized ( )

Whether NekRS itself has been initialized yet

Returns
whether NekRS is initialized

◆ isMovingMeshBoundary()

bool nekrs::isMovingMeshBoundary ( const int  boundary)

Whether the specific boundary is a moving mesh boundary

Parameters
[in]boundaryboundary ID
Returns
whether boundary is a moving mesh boundary

◆ isTemperatureBoundary()

bool nekrs::isTemperatureBoundary ( const int  boundary)

Whether the specific boundary is a specified temperature boundary

Parameters
[in]boundaryboundary ID
Returns
whether boundary is a temperature boundary

◆ limitTemperature()

void nekrs::limitTemperature ( const double *  min_T,
const double *  max_T 
)

Limit the temperature in nekRS to within the range of [min_T, max_T]

Parameters
[in]min_Tminimum temperature allowable in nekRS
[in]max_Tmaximum temperature allowable in nekRS

◆ massFlowrate()

double nekrs::massFlowrate ( const std::vector< int > &  boundary_id,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the mass flowrate over a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to compute the mass flowrate
[in]pp_meshwhich NekRS mesh to operate on
Returns
mass flowrate

◆ NboundaryFaces()

int nekrs::NboundaryFaces ( )

Total number of element faces on a boundary of the nekRS mesh summed over all processes

Returns
number of boundary element faces

◆ NboundaryID()

int nekrs::NboundaryID ( )

Number of boundary IDs in the nekRS mesh

Returns
number of boundary IDs

◆ Nelements()

int nekrs::Nelements ( )

Total number of volume elements in nekRS mesh summed over all processes

Returns
number of volume elements

◆ nestedElementsOnFace()

std::vector< std::vector< int > > nekrs::nestedElementsOnFace ( const int &  n)

Get all the element IDs that touch each face of a parent NekRS element

Parameters
[in]nNekRS polynomial order
Returns
nested element IDs, indexed by face

◆ Nfaces()

int nekrs::Nfaces ( )

Number of faces per element; because NekRS only supports HEX20, this should be 6

Returns
number of faces per mesh element

◆ NfaceVertices()

int nekrs::NfaceVertices ( )

Number of vertices required to define an element face Vertices refer to the points required to place the "corners" of an element face, and not the quadrature points. For instance, for hexahedral elements, the number of vertices per face is 4 regardless of the polynomial order.

Returns
Number of vertices per element face

◆ normalizeFlux()

bool nekrs::normalizeFlux ( const NekBoundaryCoupling nek_boundary_coupling,
const std::vector< int > &  boundary,
const double  moose_integral,
double  nek_integral,
double &  normalized_nek_integral 
)

Normalize the flux sent to nekRS to conserve the total flux

Parameters
[in]nek_boundary_couplingdata structure holding boundary coupling info
[in]boundaryboundaries for which to normalize the flux
[in]moose_integraltotal integrated flux from MOOSE to conserve
[in]nek_integraltotal integrated flux in nekRS to adjust
[out]normalized_nek_integralfinal normalized nek flux integral
Returns
whether normalization was successful, i.e. normalized_nek_integral equals moose_integral

◆ normalizeFluxBySideset()

bool nekrs::normalizeFluxBySideset ( const NekBoundaryCoupling nek_boundary_coupling,
const std::vector< int > &  boundary,
const std::vector< double > &  moose_integral,
std::vector< double > &  nek_integral,
double &  normalized_nek_integral 
)

Normalize the flux sent to nekRS to conserve the total flux

Parameters
[in]nek_boundary_couplingdata structure holding boundary coupling info
[in]boundaryboundaries for which to normalize the flux
[in]moose_integraltotal integrated flux from MOOSE to conserve
[in]nek_integraltotal integrated flux in nekRS to adjust
[out]normalized_nek_integralfinal normalized nek flux integral
Returns
whether normalization was successful, i.e. normalized_nek_integral equals moose_integral

◆ polynomialOrder()

int nekrs::polynomialOrder ( )

Polynomial order used in nekRS solution

Returns
polynomial order

◆ Pr()

double nekrs::Pr ( )

Get the Prandtl number; note that for dimensional cases, this is only guaranteed to be correct if the density, viscosity, heat capacity, and conductivity are constant.

Returns
constant Prandtl number

◆ pressure()

double nekrs::pressure ( const int  id)

Get the pressure solution at given GLL index

Parameters
[in]idGLL index
Returns
pressure value at index

◆ pressureSurfaceForce()

double nekrs::pressureSurfaceForce ( const std::vector< int > &  boundary_id,
const Point &  direction,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the integral of pressure on a surface, multiplied by the unit normal of the surface with a specified direction vector.

Parameters
[in]boundary_idNekRS boundary IDs for which to perform the integral
[in]directionunit vector to dot with the boundary surface normal
[in]pp_meshwhich NekRS mesh to operate on
Returns
pressure surface force, along a particular direction

◆ referenceArea()

double nekrs::referenceArea ( )

Get the reference area scale

Returns
reference area scale

◆ referenceFlux()

double nekrs::referenceFlux ( )

Get the reference heat flux scale, \(\rho C_pU\Delta T\)

Returns
reference heat flux scale

◆ referenceLength()

double nekrs::referenceLength ( )

Get the reference length scale

Returns
reference length scale

◆ referenceSource()

double nekrs::referenceSource ( )

Get the reference heat source scale, \(\rho C_pU\Delta T/L\)

Returns
reference heat source scale

◆ referenceVolume()

double nekrs::referenceVolume ( )

Get the reference volume scale

Returns
reference volume scale

◆ resolveType()

template<typename T >
MPI_Datatype nekrs::resolveType ( )

◆ scalar01()

double nekrs::scalar01 ( const int  id)

Get the scalar01 solution at given GLL index.

Parameters
[in]idGLL index
Returns
scalar01 value at index

◆ scalar02()

double nekrs::scalar02 ( const int  id)

Get the scalar02 solution at given GLL index.

Parameters
[in]idGLL index
Returns
scalar02 value at index

◆ scalar03()

double nekrs::scalar03 ( const int  id)

Get the scalar03 solution at given GLL index.

Parameters
[in]idGLL index
Returns
scalar03 value at index

◆ scalarFieldOffset()

int nekrs::scalarFieldOffset ( )

Offset increment for indexing into multi-volume arrays for the scalar fields. This assumes that all scalars are the same length as the temperature scalar. TODO: evaluate whether this works if nekRS uses CHT

Returns
scalar field offset

◆ scaleUsrwrk()

void nekrs::scaleUsrwrk ( const unsigned int &  slot,
const dfloat &  value 
)

Scale a slot in the usrwrk by a fixed value (multiplication)

Parameters
[in]slotslot in usrwrk to modify
[in]valuevalue to multiply on scratch slot

◆ scratchAvailable()

bool nekrs::scratchAvailable ( )

Whether the scratch space has already been allocated by the user

Returns
whether scratch space is already allocated

◆ setAbsoluteTol()

void nekrs::setAbsoluteTol ( double  tol)

Set the absolute tolerance for checking energy conservation in data transfers to Nek

Parameters
[in]toltolerance

◆ setNekSetupTime()

void nekrs::setNekSetupTime ( const double &  time)

Nek's runtime statistics are formed by collecting a timer of both the initialization and accumulated run time. We unfortunately have to split this across multiple classes, so if we want correct times we need to have NekInitAction save the value of the time spent on initialization.

Parameters
[in]timetime spent on initialization

◆ setRelativeTol()

void nekrs::setRelativeTol ( double  tol)

Set the relative tolerance for checking energy conservation in data transfers to Nek

Parameters
[in]toltolerance

◆ setStartTime()

void nekrs::setStartTime ( const double &  start)

Set the start time used by NekRS

Parameters
[in]startstart time

◆ sideExtremeValue()

double nekrs::sideExtremeValue ( const std::vector< int > &  boundary_id,
const field::NekFieldEnum field,
const nek_mesh::NekMeshEnum  pp_mesh,
const bool  max 
)

Find the extreme of a given field over a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to find the extreme value
[in]fieldfield to find the maximum value of
[in]pp_meshwhich NekRS mesh to operate on
[in]maxwhether to take the maximum (or if false, the minimum)
Returns
max or min value of field on boundary

◆ sideIntegral()

double nekrs::sideIntegral ( const std::vector< int > &  boundary_id,
const field::NekFieldEnum integrand,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the area integral of a given integrand over a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to perform the integral
[in]integrandfield to integrate
[in]pp_meshwhich NekRS mesh to operate on
Returns
area integral of a field

◆ sideMassFluxWeightedIntegral()

double nekrs::sideMassFluxWeightedIntegral ( const std::vector< int > &  boundary_id,
const field::NekFieldEnum integrand,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the mass flux weighted integral of a given integrand over a set of boundary IDs

Parameters
[in]boundary_idnekRS boundary IDs for which to perform the integral
[in]integrandfield to integrate and weight by mass flux
[in]pp_meshwhich NekRS mesh to operate on
Returns
mass flux weighted area average of a field

◆ storeBoundaryCoupling()

void nekrs::storeBoundaryCoupling ( const std::vector< int > &  boundary_id,
int &  N 
)

Store the rank-local element, element-local face, and rank ownership for boundary coupling

Parameters
[in]boundary_idboundaries through which nekRS will be coupled
[out]Ntotal number of surface elements

◆ temperature()

double nekrs::temperature ( const int  id)

Get the temperature solution at given GLL index.

Because nekRS stores all the passive scalars together in one flat array, this routine simply indices into the entire passive scalar solution. In order to get temperature, you should only index up to nrs->cds->fieldOffset.

Parameters
[in]idGLL index
Returns
temperature value at index

◆ temperatureBoundaryType()

const std::string nekrs::temperatureBoundaryType ( const int  boundary)

String name indicating the temperature boundary condition type on a given boundary

Parameters
[in]boundaryboundary ID
Returns
string name of boundary condition type

◆ temperatureMesh()

mesh_t* nekrs::temperatureMesh ( )

Get the mesh for the temperature scalar

Returns
temperature mesh

◆ unity()

double nekrs::unity ( const int  id)

Return unity, for cases where the integrand or operator we are generalizing acts on 1

Parameters
[in]idGLL index
Returns
unity

◆ updateHostMeshParameters()

void nekrs::updateHostMeshParameters ( )

Update the mesh parameters on host.

◆ usrwrkSideIntegral()

std::vector<double> nekrs::usrwrkSideIntegral ( const unsigned int &  slot,
const std::vector< int > &  boundary,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Integrate the scratch space over boundaries

Parameters
[in]slotslot in scratch space
[in]boundaryboundaries over which to integrate the scratch space
[in]pp_meshportion of NekRS mesh to integrate over
Returns
boundary integrated scratch space, with one value per sideset

◆ usrwrkVolumeIntegral()

double nekrs::usrwrkVolumeIntegral ( const unsigned int &  slot,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Volume integrate the scratch space

Parameters
[in]slotslot in scratch space to i ntegrat
[in]pp_meshNekRS mesh to integrate over
Returns
volume integrated scratch space

◆ validBoundaryIDs()

bool nekrs::validBoundaryIDs ( const std::vector< int > &  boundary_id,
int &  first_invalid_id,
int &  n_boundaries 
)

Whether the provided boundary IDs are all valid in the nekRS mesh

Parameters
[in]boundary_idvector of boundary IDs to check
[out]first_invalid_idfirst invalid ID encountered for printing an error on the MOOSE side
[out]n_boundariesmaximum valid boundary ID for printing an error on the MOOSE side
Returns
whether all boundaries are valid

◆ velocity()

double nekrs::velocity ( const int  id)

Get the magnitude of the velocity solution at given GLL index

Parameters
[in]idGLL index
Returns
velocity magnitude at index

◆ velocity_x()

double nekrs::velocity_x ( const int  id)

Get the x-velocity at given GLL index

Parameters
[in]idGLL index
Returns
x-velocity at index

◆ velocity_y()

double nekrs::velocity_y ( const int  id)

Get the y-velocity at given GLL index

Parameters
[in]idGLL index
Returns
y-velocity at index

◆ velocity_z()

double nekrs::velocity_z ( const int  id)

Get the z-velocity at given GLL index

Parameters
[in]idGLL index
Returns
z-velocity at index

◆ velocityFieldOffset()

int nekrs::velocityFieldOffset ( )

Offset increment for indexing into the velocity array

Returns
velocity field offset

◆ viscosity()

double nekrs::viscosity ( )

Get the viscosity used in the definition of the Reynolds number; note that for dimensional cases, this is only guaranteed to be correct if the viscosity is constant.

Returns
constant dynamic viscosity

◆ volume()

double nekrs::volume ( const nek_mesh::NekMeshEnum  pp_mesh)

Compute the volume over the entire scalar mesh

Parameters
[in]pp_meshwhich NekRS mesh to operate on
Returns
volume integral

◆ volumeExtremeValue()

double nekrs::volumeExtremeValue ( const field::NekFieldEnum field,
const nek_mesh::NekMeshEnum  pp_mesh,
const bool  max 
)

Find the extreme value of a given field over the entire nekRS domain

Parameters
[in]fieldfield to find the minimum value of
[in]pp_meshwhich NekRS mesh to operate on
[in]maxwhether to take the maximum (or if false, the minimum)
Returns
max or min value of field in volume

◆ volumeIntegral()

double nekrs::volumeIntegral ( const field::NekFieldEnum integrand,
const double &  volume,
const nek_mesh::NekMeshEnum  pp_mesh 
)

Compute the volume integral of a given integrand over the entire scalar mesh

Parameters
[in]integrandfield to integrate
[in]volumevolume of the domain (only used for dimensionalizing temperature)
[in]pp_meshwhich NekRS mesh to operate on
Returns
volume integral of a field

◆ write_field_file()

void nekrs::write_field_file ( const std::string &  prefix,
const dfloat  time,
const int &  step 
)

Write a field file containing pressure, velocity, and scalars with given prefix

Parameters
[in]prefixthree-character prefix
[in]timetime
[in]steptime step index

◆ write_usrwrk_field_file()

void nekrs::write_usrwrk_field_file ( const int &  slot,
const std::string &  prefix,
const dfloat &  time,
const int &  step,
const bool &  write_coords 
)

Write a field file containing a specific slot of the nrs->usrwrk scratch space; this will write the field to the 'temperature' slot in a field file.

Parameters
[in]slotindex in the nrs->usrwrk array to write
[in]prefixprefix for file name
[in]timesimulation time to write file for
[in]steptime step index
[in]write_coordswhether to write the mesh coordinates

◆ x_displacement()

void nekrs::x_displacement ( const int  id,
const dfloat  value 
)

Write a value into the x-displacement

Parameters
[in]idindex
[in]valuevalue to write

◆ y_displacement()

void nekrs::y_displacement ( const int  id,
const dfloat  value 
)

Write a value into the y-displacement

Parameters
[in]idindex
[in]valuevalue to write

◆ z_displacement()

void nekrs::z_displacement ( const int  id,
const dfloat  value 
)

Write a value into the z-displacement

Parameters
[in]idindex
[in]valuevalue to write

Variable Documentation

◆ solutionPointer

void(*)(int, dfloat) nekrs::solutionPointer(const field::NekWriteEnum &field)

Get pointer to various solution functions (for reading only) based on enumeration

Parameters
[in]fieldfield to return a pointer to
Returns
function pointer to method that returns said field as a function of GLL index

Write various solution functions based on enumeration

Parameters
[in]fieldfield to write