Cardinal
|
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) |
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.
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
[in] | base_counts | once multiplied by 'multiplier', the number of counts on each rank |
[in] | input | rank-local data |
[out] | output | collected result |
[in] | multiplier | constant multiplier to set on each count indicator |
double nekrs::area | ( | const std::vector< int > & | boundary_id, |
const nek_mesh::NekMeshEnum | pp_mesh | ||
) |
Compute the area of a set of boundary IDs
[in] | boundary_id | nekRS boundary IDs for which to perform the integral |
[in] | pp_mesh | which NekRS mesh to operate on |
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).
[in] | buildOnly | whether NekRS is to be run in build-only mode |
int nekrs::buildOnly | ( | ) |
Whether NekRS was run in JIT build-only mode
Point nekrs::centroid | ( | int | local_elem_id | ) |
Compute the centroid given a local element ID (NOTE: returns in dimensional form)
[in] | local_elem_id | local element ID on this rank |
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)
[in] | local_elem_id | local element ID on this rank |
[in] | local_face_id | local face ID on the element |
int nekrs::commRank | ( | ) |
Get the process rank
int nekrs::commSize | ( | ) |
Get the communicator size
void nekrs::copyDeformationToDevice | ( | ) |
Copy the deformation from host to device.
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:
[in] | n | order of mesh (1 = first-order, 2 = second-order) |
[in] | exact | whether the MOOSE elements will exactly represent the NekRS mesh |
int nekrs::dim | ( | ) |
Mesh dimension
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.
[in] | field | physical interpretation of value to dimensionalize |
[out] | value | value to dimensionalize |
void nekrs::dimensionalizeArea | ( | double & | integral | ) |
Dimensionalize an area
[in] | integral | integral to dimensionalize |
void nekrs::dimensionalizeSideIntegral | ( | const field::NekFieldEnum & | integrand, |
const Real & | area, | ||
double & | integral | ||
) |
Dimensionalize a given integral of f over a side, i.e. fdS
[in] | integrand | field to dimensionalize |
[in] | area | area of the boundary |
[in] | integral | integral to dimensionalize |
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
[in] | integrand | field to dimensionalize |
[in] | boundary_id | boundary IDs for the integral |
[in] | integral | integral to dimensionalize |
[in] | pp_mesh | which NekRS mesh to operate on |
void nekrs::dimensionalizeVolume | ( | double & | integral | ) |
Dimensionalize a volume
[in] | integral | integral to dimensionalize |
void nekrs::dimensionalizeVolumeIntegral | ( | const field::NekFieldEnum & | integrand, |
const Real & | volume, | ||
double & | integral | ||
) |
Dimensionalize a given integral of f over volume, i.e. fdV
[in] | integrand | field to dimensionalize |
[in] | volume | volume of the domain (only used for dimensionalizing temperature) |
[in] | integral | integral to dimensionalize |
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
[in] | base_counts | unit-wise receiving counts for each process |
[out] | counts | receiving counts from each process |
[out] | displacement | displacement for each process's counts |
[in] | multiplier | optional multiplier on the face-based data |
bool nekrs::endControlElapsedTime | ( | ) |
Whether nekRS's input file intends to terminate the simulation based on a wall time
bool nekrs::endControlNumSteps | ( | ) |
Whether nekRS's input file intends to terminate the simulation based on a number of steps
bool nekrs::endControlTime | ( | ) |
Whether nekRS's input file intends to terminate the simulation based on an end time
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.
mesh_t* nekrs::flowMesh | ( | ) |
Get the mesh for the flow solve
void nekrs::flux | ( | const int | id, |
const dfloat | value | ||
) |
Write a value into the user scratch space that holds the flux
[in] | id | index |
[in] | value | value to write |
void nekrs::freeScratch | ( | ) |
Free the scratch space.
mesh_t* nekrs::getMesh | ( | const nek_mesh::NekMeshEnum | pp_mesh | ) |
Get the mesh to act on
[in] | pp_mesh | which NekRS mesh to operate on |
double nekrs::getNekSetupTime | ( | ) |
Get time spent on initialization
dfloat* nekrs::getSgeo | ( | ) |
dfloat* nekrs::getVgeo | ( | ) |
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)
[in] | local_elem_id | local element ID on this rank |
[in] | local_node_id | local node ID on this element |
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)
[in] | local_elem_id | local element ID on this rank |
[in] | local_face_id | local face ID on this element |
[in] | local_node_id | local node ID on this element |
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
[in] | offset | in the gradient field for each component (grad_x, grad_y, or grad_z) |
[in] | f | field to compute the gradient of |
[in] | pp_mesh | which NekRS mesh to operate on |
[out] | grad_f | gradient of field |
bool nekrs::hasBlendingSolver | ( | ) |
Whether nekRS's input file has the blending mesh solver
bool nekrs::hasCHT | ( | ) |
Whether nekRS's input file has CHT
bool nekrs::hasHeatSourceKernel | ( | ) |
Whether nekRS contains an OCCA kernel to apply a source to the passive scalar equations
bool nekrs::hasMovingMesh | ( | ) |
Whether nekRS's input file indicates a moving mesh
bool nekrs::hasScalarVariable | ( | int | scalarId | ) |
Whether nekRS's input file indicates that the problem has a scalar0(scalarId) variable
[in] | scalarId | scalar number, i.e. for scalar03 scalarId=3 |
bool nekrs::hasTemperatureSolve | ( | ) |
Whether nekRS actually solves for temperature (as opposed to setting its solver to 'none')
bool nekrs::hasTemperatureVariable | ( | ) |
Whether nekRS's input file indicates that the problem has a temperature variable
bool nekrs::hasUserMeshSolver | ( | ) |
Whether nekRS's input file has the user mesh solver
bool nekrs::hasVariableDt | ( | ) |
Whether nekRS's input file indicates a variable time stepping scheme
void nekrs::heat_source | ( | const int | id, |
const dfloat | value | ||
) |
Write a value into the user scratch space that holds the volumetric heat source
[in] | id | index |
[in] | value | value to write |
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
[in] | boundary_id | nekRS boundary IDs for which to perform the integral |
[in] | pp_mesh | which NekRS mesh to operate on |
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
[in] | U_ref | reference velocity |
[in] | T_ref | reference temperature |
[in] | dT_ref | reference temperature range |
[in] | L_ref | reference length scale |
[in] | rho_ref | reference density |
[in] | Cp_ref | reference heat capacity |
void nekrs::initializeHostMeshParameters | ( | ) |
Allocate memory for the host mesh parameters.
void nekrs::initializeScratch | ( | const unsigned int & | n_slots | ) |
Initialize scratch space for data to get sent into NekRS
[in] | n_slots | number of slots (for volume arrays) to allocate |
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
[in] | scratch | available scratch space for the calculation |
[in] | I | interpolation matrix |
[in] | x | face data to be interpolated |
[in] | N | number of points in 1-D to be interpolated |
[out] | Ix | interpolated data |
[in] | M | resulting number of interpolated points in 1-D |
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
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.
[out] | I | interpolation matrix |
[in] | starting_points | number of points in the source quadrature rule |
[in] | ending_points | number of points in the end quadrature rule |
bool nekrs::isHeatFluxBoundary | ( | const int | boundary | ) |
Whether the specific boundary is a flux boundary
[in] | boundary | boundary ID |
bool nekrs::isInitialized | ( | ) |
Whether NekRS itself has been initialized yet
bool nekrs::isMovingMeshBoundary | ( | const int | boundary | ) |
Whether the specific boundary is a moving mesh boundary
[in] | boundary | boundary ID |
bool nekrs::isTemperatureBoundary | ( | const int | boundary | ) |
Whether the specific boundary is a specified temperature boundary
[in] | boundary | boundary ID |
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]
[in] | min_T | minimum temperature allowable in nekRS |
[in] | max_T | maximum temperature allowable in nekRS |
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
[in] | boundary_id | nekRS boundary IDs for which to compute the mass flowrate |
[in] | pp_mesh | which NekRS mesh to operate on |
int nekrs::NboundaryFaces | ( | ) |
Total number of element faces on a boundary of the nekRS mesh summed over all processes
int nekrs::NboundaryID | ( | ) |
Number of boundary IDs in the nekRS mesh
int nekrs::Nelements | ( | ) |
Total number of volume elements in nekRS mesh summed over all processes
std::vector< std::vector< int > > nekrs::nestedElementsOnFace | ( | const int & | n | ) |
Get all the element IDs that touch each face of a parent NekRS element
[in] | n | NekRS polynomial order |
int nekrs::Nfaces | ( | ) |
Number of faces per element; because NekRS only supports HEX20, this should be 6
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.
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
[in] | nek_boundary_coupling | data structure holding boundary coupling info |
[in] | boundary | boundaries for which to normalize the flux |
[in] | moose_integral | total integrated flux from MOOSE to conserve |
[in] | nek_integral | total integrated flux in nekRS to adjust |
[out] | normalized_nek_integral | final normalized nek flux integral |
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
[in] | nek_boundary_coupling | data structure holding boundary coupling info |
[in] | boundary | boundaries for which to normalize the flux |
[in] | moose_integral | total integrated flux from MOOSE to conserve |
[in] | nek_integral | total integrated flux in nekRS to adjust |
[out] | normalized_nek_integral | final normalized nek flux integral |
int nekrs::polynomialOrder | ( | ) |
Polynomial order used in nekRS solution
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.
double nekrs::pressure | ( | const int | id | ) |
Get the pressure solution at given GLL index
[in] | id | GLL index |
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. This represents the force that the fluid exerts ON the boundary.
[in] | boundary_id | NekRS boundary IDs for which to perform the integral |
[in] | direction | unit vector to dot with the boundary surface normal |
[in] | pp_mesh | which NekRS mesh to operate on |
double nekrs::referenceArea | ( | ) |
Get the reference area scale
double nekrs::referenceFlux | ( | ) |
Get the reference heat flux scale, \(\rho C_pU\Delta T\)
double nekrs::referenceLength | ( | ) |
Get the reference length scale
double nekrs::referenceSource | ( | ) |
Get the reference heat source scale, \(\rho C_pU\Delta T/L\)
double nekrs::referenceVolume | ( | ) |
Get the reference volume scale
MPI_Datatype nekrs::resolveType | ( | ) |
double nekrs::scalar01 | ( | const int | id | ) |
Get the scalar01 solution at given GLL index.
[in] | id | GLL index |
double nekrs::scalar02 | ( | const int | id | ) |
Get the scalar02 solution at given GLL index.
[in] | id | GLL index |
double nekrs::scalar03 | ( | const int | id | ) |
Get the scalar03 solution at given GLL index.
[in] | id | GLL index |
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
void nekrs::scaleUsrwrk | ( | const unsigned int & | slot, |
const dfloat & | value | ||
) |
Scale a slot in the usrwrk by a fixed value (multiplication)
[in] | slot | slot in usrwrk to modify |
[in] | value | value to multiply on scratch slot |
bool nekrs::scratchAvailable | ( | ) |
Whether the scratch space has already been allocated by the user
void nekrs::setAbsoluteTol | ( | double | tol | ) |
Set the absolute tolerance for checking energy conservation in data transfers to Nek
[in] | tol | tolerance |
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.
[in] | time | time spent on initialization |
void nekrs::setRelativeTol | ( | double | tol | ) |
Set the relative tolerance for checking energy conservation in data transfers to Nek
[in] | tol | tolerance |
void nekrs::setStartTime | ( | const double & | start | ) |
Set the start time used by NekRS
[in] | start | start time |
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
[in] | boundary_id | nekRS boundary IDs for which to find the extreme value |
[in] | field | field to find the maximum value of |
[in] | pp_mesh | which NekRS mesh to operate on |
[in] | max | whether to take the maximum (or if false, the minimum) |
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
[in] | boundary_id | nekRS boundary IDs for which to perform the integral |
[in] | integrand | field to integrate |
[in] | pp_mesh | which NekRS mesh to operate on |
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
[in] | boundary_id | nekRS boundary IDs for which to perform the integral |
[in] | integrand | field to integrate and weight by mass flux |
[in] | pp_mesh | which NekRS mesh to operate on |
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
[in] | boundary_id | boundaries through which nekRS will be coupled |
[out] | N | total number of surface elements |
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.
[in] | id | GLL index |
const std::string nekrs::temperatureBoundaryType | ( | const int | boundary | ) |
String name indicating the temperature boundary condition type on a given boundary
[in] | boundary | boundary ID |
mesh_t* nekrs::temperatureMesh | ( | ) |
Get the mesh for the temperature scalar
double nekrs::unity | ( | const int | id | ) |
Return unity, for cases where the integrand or operator we are generalizing acts on 1
[in] | id | GLL index |
void nekrs::updateHostMeshParameters | ( | ) |
Update the mesh parameters on host.
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
[in] | slot | slot in scratch space |
[in] | boundary | boundaries over which to integrate the scratch space |
[in] | pp_mesh | portion of NekRS mesh to integrate over |
double nekrs::usrwrkVolumeIntegral | ( | const unsigned int & | slot, |
const nek_mesh::NekMeshEnum | pp_mesh | ||
) |
Volume integrate the scratch space
[in] | slot | slot in scratch space to i ntegrat |
[in] | pp_mesh | NekRS mesh to integrate over |
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
[in] | boundary_id | vector of boundary IDs to check |
[out] | first_invalid_id | first invalid ID encountered for printing an error on the MOOSE side |
[out] | n_boundaries | maximum valid boundary ID for printing an error on the MOOSE side |
double nekrs::velocity | ( | const int | id | ) |
Get the magnitude of the velocity solution at given GLL index
[in] | id | GLL index |
double nekrs::velocity_x | ( | const int | id | ) |
Get the x-velocity at given GLL index
[in] | id | GLL index |
double nekrs::velocity_y | ( | const int | id | ) |
Get the y-velocity at given GLL index
[in] | id | GLL index |
double nekrs::velocity_z | ( | const int | id | ) |
Get the z-velocity at given GLL index
[in] | id | GLL index |
int nekrs::velocityFieldOffset | ( | ) |
Offset increment for indexing into the velocity array
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.
double nekrs::volume | ( | const nek_mesh::NekMeshEnum | pp_mesh | ) |
Compute the volume over the entire scalar mesh
[in] | pp_mesh | which NekRS mesh to operate on |
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
[in] | field | field to find the minimum value of |
[in] | pp_mesh | which NekRS mesh to operate on |
[in] | max | whether to take the maximum (or if false, the minimum) |
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
[in] | integrand | field to integrate |
[in] | volume | volume of the domain (only used for dimensionalizing temperature) |
[in] | pp_mesh | which NekRS mesh to operate on |
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
[in] | prefix | three-character prefix |
[in] | time | time |
[in] | step | time step index |
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.
[in] | slot | index in the nrs->usrwrk array to write |
[in] | prefix | prefix for file name |
[in] | time | simulation time to write file for |
[in] | step | time step index |
[in] | write_coords | whether to write the mesh coordinates |
void nekrs::x_displacement | ( | const int | id, |
const dfloat | value | ||
) |
Write a value into the x-displacement
[in] | id | index |
[in] | value | value to write |
void nekrs::y_displacement | ( | const int | id, |
const dfloat | value | ||
) |
Write a value into the y-displacement
[in] | id | index |
[in] | value | value to write |
void nekrs::z_displacement | ( | const int | id, |
const dfloat | value | ||
) |
Write a value into the z-displacement
[in] | id | index |
[in] | value | value to write |
void(*)(int, dfloat) nekrs::solutionPointer(const field::NekWriteEnum &field) |
Get pointer to various solution functions (for reading only) based on enumeration
[in] | field | field to return a pointer to |
Write various solution functions based on enumeration
[in] | field | field to write |