Go to the documentation of this file.
22 #include "MooseTypes.h"
33 #include "libmesh/point.h"
48 static int build_only;
107 void write_usrwrk_field_file(
const int & slot,
const std::string & prefix,
const dfloat & time,
const int & step,
const bool & write_coords);
115 void write_field_file(
const std::string & prefix,
const dfloat time,
const int & step);
294 template <
typename T>
295 void allgatherv(
const std::vector<int> & base_counts,
298 const int multiplier = 1);
310 const int multiplier);
331 double * scratch,
const double * I,
double * x,
int N,
double * Ix,
int M);
340 Point
centroidFace(
int local_elem_id,
int local_face_id);
355 Point
gllPoint(
int local_elem_id,
int local_node_id);
365 Point
gllPointFace(
int local_elem_id,
int local_face_id,
int local_node_id);
375 const std::vector<int> & boundary,
391 void scaleUsrwrk(
const unsigned int & slot,
const dfloat & value);
403 const std::vector<int> & boundary,
404 const std::vector<double> & moose_integral,
405 std::vector<double> & nek_integral,
406 double & normalized_nek_integral);
418 const std::vector<int> & boundary,
419 const double moose_integral,
421 double & normalized_nek_integral);
488 const std::vector<int> & boundary_id,
509 double massFlowrate(
const std::vector<int> & boundary_id,
556 void gradient(
const int offset,
const double * f,
double * grad_f,
662 validBoundaryIDs(
const std::vector<int> & boundary_id,
int & first_invalid_id,
int & n_boundaries);
794 double unity(
const int id);
829 void flux(
const int id,
const dfloat value);
836 void heat_source(
const int id,
const dfloat value);
872 const double rho_ref,
873 const double Cp_ref);
921 template <
typename T>
931 template <
typename T>
933 allgatherv(
const std::vector<int> & base_counts,
const T * input, T * output,
const int multiplier)
935 int * recvCounts = (
int *)calloc(
commSize(),
sizeof(int));
936 int * displacement = (
int *)calloc(
commSize(),
sizeof(int));
939 MPI_Allgatherv(input,
943 (
const int *)recvCounts,
944 (
const int *)displacement,
946 platform->comm.mpiComm);
int mesh_velocity_y
y-velocity of moving boundary (for mesh blending solver)
Definition: NekInterface.h:687
bool endControlNumSteps()
double(*)(int) solutionPointer(const field::NekFieldEnum &field)
Definition: NekInterface.h:739
double source_ref
Definition: NekInterface.h:729
bool isHeatFluxBoundary(const int boundary)
Point centroidFace(int local_elem_id, int local_face_id)
void interpolationMatrix(double *I, int starting_points, int ending_points)
double scalar01(const int id)
Get the scalar01 solution at given GLL index.
Point centroid(int local_elem_id)
MPI_Datatype resolveType()
NekWriteEnum
Enumeration of possible fields to write in nekRS.
Definition: CardinalEnums.h:104
void storeBoundaryCoupling(const std::vector< int > &boundary_id, int &N)
Definition: NekBoundaryCoupling.h:29
std::vector< double > usrwrkSideIntegral(const unsigned int &slot, const std::vector< int > &boundary, const nek_mesh::NekMeshEnum pp_mesh)
void flux(const int id, const dfloat value)
void displacementAndCounts(const std::vector< int > &base_counts, int *counts, int *displacement, const int multiplier)
Definition: CardinalEnums.h:85
void heat_source(const int id, const dfloat value)
double velocity(const int id)
bool hasHeatSourceKernel()
double dT_ref
Definition: NekInterface.h:715
bool validBoundaryIDs(const std::vector< int > &boundary_id, int &first_invalid_id, int &n_boundaries)
double sideExtremeValue(const std::vector< int > &boundary_id, const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
double massFlowrate(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
double usrwrkVolumeIntegral(const unsigned int &slot, const nek_mesh::NekMeshEnum pp_mesh)
bool nondimensional_T
Definition: NekInterface.h:731
int heat_source
volumetric heat source (for volumetric heating)
Definition: NekInterface.h:681
double Cp_ref
Definition: NekInterface.h:725
@ max
Definition: CardinalEnums.h:122
void interpolateSurfaceFaceHex3D(double *scratch, const double *I, double *x, int N, double *Ix, int M)
double area(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
void setAbsoluteTol(double tol)
double T_ref
Definition: NekInterface.h:713
double pressureSurfaceForce(const std::vector< int > &boundary_id, const Point &direction, const nek_mesh::NekMeshEnum pp_mesh)
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)
double scalar03(const int id)
Get the scalar03 solution at given GLL index.
double A_ref
Definition: NekInterface.h:719
bool hasTemperatureVariable()
double volume(const nek_mesh::NekMeshEnum pp_mesh)
double U_ref
Definition: NekInterface.h:711
Characteristic scales assumed in nekRS if using a non-dimensional solution.
Definition: NekInterface.h:709
Point gllPointFace(int local_elem_id, int local_face_id, int local_node_id)
int boundary_scalar02
boundary scalar02 (for separate domain coupling)
Definition: NekInterface.h:702
mesh_t * getMesh(const nek_mesh::NekMeshEnum pp_mesh)
void y_displacement(const int id, const dfloat value)
void setRelativeTol(double tol)
double pressure(const int id)
void updateHostMeshParameters()
Update the mesh parameters on host.
void dimensionalizeSideIntegral(const field::NekFieldEnum &integrand, const Real &area, double &integral)
void allgatherv(const std::vector< int > &base_counts, const T *input, T *output, const int multiplier=1)
Definition: NekInterface.h:933
double volumeExtremeValue(const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
void write_field_file(const std::string &prefix, const dfloat time, const int &step)
double velocity_x(const int id)
int velocityFieldOffset()
void limitTemperature(const double *min_T, const double *max_T)
void initializeScratch(const unsigned int &n_slots)
double rho_ref
Definition: NekInterface.h:723
NekMeshEnum
Definition: CardinalEnums.h:44
double L_ref
Definition: NekInterface.h:717
Point gllPoint(int local_elem_id, int local_node_id)
void copyDeformationToDevice()
Copy the deformation from host to device.
void buildOnly(int buildOnly)
void dimensionalizeVolumeIntegral(const field::NekFieldEnum &integrand, const Real &volume, double &integral)
int NfaceVertices()
Number of vertices required to define an element face Vertices refer to the points required to place ...
bool hasScalarVariable(int scalarId)
void interpolateVolumeHex3D(const double *I, double *x, int N, double *Ix, int M)
double heatFluxIntegral(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
int flux
boundary heat flux (for conjugate heat transfer)
Definition: NekInterface.h:678
bool hasTemperatureSolve()
void dimensionalize(const field::NekFieldEnum &field, double &value)
Dimensionalize a field by multiplying the nondimensional form by the reference.
double volumeIntegral(const field::NekFieldEnum &integrand, const double &volume, const nek_mesh::NekMeshEnum pp_mesh)
int boundary_velocity
boundary velocity (for separate domain coupling)
Definition: NekInterface.h:693
double unity(const int id)
double velocity_z(const int id)
double sideMassFluxWeightedIntegral(const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
Definition: NekInterface.h:675
double temperature(const int id)
Get the temperature solution at given GLL index.
void x_displacement(const int id, const dfloat value)
bool endControlElapsedTime()
double sideIntegral(const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
int boundary_scalar01
boundary scalar01 (for separate domain coupling)
Definition: NekInterface.h:699
int boundary_scalar03
boundary scalar03 (for separate domain coupling)
Definition: NekInterface.h:705
void gradient(const int offset, const double *f, double *grad_f, const nek_mesh::NekMeshEnum pp_mesh)
mesh_t * temperatureMesh()
void initializeHostMeshParameters()
Allocate memory for the host mesh parameters.
int mesh_velocity_x
x-velocity of moving boundary (for mesh blending solver)
Definition: NekInterface.h:684
double velocity_y(const int id)
void setNekSetupTime(const double &time)
void write_usrwrk_field_file(const int &slot, const std::string &prefix, const dfloat &time, const int &step, const bool &write_coords)
void dimensionalizeArea(double &integral)
double flux_ref
Definition: NekInterface.h:727
void freeScratch()
Free the scratch space.
const std::string temperatureBoundaryType(const int boundary)
bool normalizeFlux(const NekBoundaryCoupling &nek_boundary_coupling, const std::vector< int > &boundary, const double moose_integral, double nek_integral, double &normalized_nek_integral)
int boundary_temperature
boundary temperature (for separate domain coupling)
Definition: NekInterface.h:696
Cardinal-specific nekRS API.
Definition: NekUility.C:21
bool isTemperatureBoundary(const int boundary)
void setStartTime(const double &start)
double V_ref
Definition: NekInterface.h:721
bool isMovingMeshBoundary(const int boundary)
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)
void scaleUsrwrk(const unsigned int &slot, const dfloat &value)
void dimensionalizeVolume(double &integral)
NekFieldEnum
Enumeration of possible fields to read from nekRS.
Definition: CardinalEnums.h:88
int mesh_velocity_z
z-velocity of moving boundary (for mesh blending solver)
Definition: NekInterface.h:690
void z_displacement(const int id, const dfloat value)
double scalar02(const int id)
Get the scalar02 solution at given GLL index.