Loading [MathJax]/extensions/tex2jax.js
Cardinal
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
NekInterface.h
Go to the documentation of this file.
1/********************************************************************/
2/* SOFTWARE COPYRIGHT NOTIFICATION */
3/* Cardinal */
4/* */
5/* (c) 2021 UChicago Argonne, LLC */
6/* ALL RIGHTS RESERVED */
7/* */
8/* Prepared by UChicago Argonne, LLC */
9/* Under Contract No. DE-AC02-06CH11357 */
10/* With the U. S. Department of Energy */
11/* */
12/* Prepared by Battelle Energy Alliance, LLC */
13/* Under Contract No. DE-AC07-05ID14517 */
14/* With the U. S. Department of Energy */
15/* */
16/* See LICENSE for full restrictions */
17/********************************************************************/
18
19#pragma once
20
21#include "CardinalEnums.h"
22#include "MooseTypes.h"
23#include "NekBoundaryCoupling.h"
24#include "NekVolumeCoupling.h"
25
26#include "inipp.hpp"
27#include "nekrs.hpp"
28#include "bcMap.hpp"
29#include "udf.hpp"
30#include "inipp.hpp"
31#include "mesh.h"
32
33#include "libmesh/point.h"
34
35#include <string>
36#include <vector>
37
45namespace nekrs
46{
47
48static int build_only;
49
52
55
56dfloat * getSgeo();
57dfloat * getVgeo();
58
65
70void setAbsoluteTol(double tol);
71
78Real scratchUnits(const int slot);
79
84void nondimensional(const bool n);
85
90void setRelativeTol(double tol);
91
98
106void setNekSetupTime(const double & time);
107
113
118void setStartTime(const double & start);
119
125
135void write_usrwrk_field_file(const int & slot, const std::string & prefix, const dfloat & time, const int & step, const bool & write_coords);
136
143void write_field_file(const std::string & prefix, const dfloat time, const int & step);
144
152
158
162void interpolateVolumeHex3D(const double * I, double * x, int N, double * Ix, int M);
163
168bool hasCHT();
169
175
181
187
193
199
205
211
219
225
231
238mesh_t * entireMesh();
239
244mesh_t * flowMesh();
245
251
257mesh_t * getMesh(const nek_mesh::NekMeshEnum pp_mesh);
258
264
270
276
282
288bool hasScalarVariable(int scalarId);
289
295
301
306void initializeScratch(const unsigned int & n_slots);
307
310
316double viscosity();
317
323double Pr();
324
327
328template <typename T>
329void allgatherv(const std::vector<int> & base_counts,
330 const T * input,
331 T * output,
332 const int multiplier = 1);
333
341void displacementAndCounts(const std::vector<int> & base_counts,
342 int * counts,
343 int * displacement,
344 const int multiplier);
345
353void interpolationMatrix(double * I, int starting_points, int ending_points);
354
365 double * scratch, const double * I, double * x, int N, double * Ix, int M);
366
374Point centroidFace(int local_elem_id, int local_face_id);
375
381Point centroid(int local_elem_id);
382
389Point gllPoint(int local_elem_id, int local_node_id);
390
399Point gllPointFace(int local_elem_id, int local_face_id, int local_node_id);
400
408std::vector<double> usrwrkSideIntegral(const unsigned int & slot,
409 const std::vector<int> & boundary,
410 const nek_mesh::NekMeshEnum pp_mesh);
411
418double usrwrkVolumeIntegral(const unsigned int & slot, const nek_mesh::NekMeshEnum pp_mesh);
419
425void scaleUsrwrk(const unsigned int & slot, const dfloat & value);
426
436bool normalizeFluxBySideset(const NekBoundaryCoupling & nek_boundary_coupling,
437 const std::vector<int> & boundary,
438 const std::vector<double> & moose_integral,
439 std::vector<double> & nek_integral,
440 double & normalized_nek_integral);
441
451bool normalizeFlux(const NekBoundaryCoupling & nek_boundary_coupling,
452 const std::vector<int> & boundary,
453 const double moose_integral,
454 double nek_integral,
455 double & normalized_nek_integral);
456
463double area(const std::vector<int> & boundary_id, const nek_mesh::NekMeshEnum pp_mesh);
464
472double sideIntegral(const std::vector<int> & boundary_id, const field::NekFieldEnum & integrand,
473 const nek_mesh::NekMeshEnum pp_mesh);
474
480double volume(const nek_mesh::NekMeshEnum pp_mesh);
481
486void dimensionalizeVolume(double & integral);
487
492void dimensionalizeArea(double & integral);
493
501 const Real & volume,
502 double & integral);
503
511 const Real & area,
512 double & integral);
513
522 const std::vector<int> & boundary_id,
523 double & integral,
524 const nek_mesh::NekMeshEnum pp_mesh);
525
533double volumeIntegral(const field::NekFieldEnum & integrand,
534 const double & volume,
535 const nek_mesh::NekMeshEnum pp_mesh);
536
543double massFlowrate(const std::vector<int> & boundary_id,
544 const nek_mesh::NekMeshEnum pp_mesh);
545
553double sideMassFluxWeightedIntegral(const std::vector<int> & boundary_id,
554 const field::NekFieldEnum & integrand,
555 const nek_mesh::NekMeshEnum pp_mesh);
556
566double pressureSurfaceForce(const std::vector<int> & boundary_id, const Point & direction, const nek_mesh::NekMeshEnum pp_mesh);
567
574double heatFluxIntegral(const std::vector<int> & boundary_id,
575 const nek_mesh::NekMeshEnum pp_mesh);
576
582void limitTemperature(const double * min_T, const double * max_T);
583
591void gradient(const int offset, const double * f, double * grad_f,
592 const nek_mesh::NekMeshEnum pp_mesh);
593
602 const nek_mesh::NekMeshEnum pp_mesh,
603 const bool max);
604
613double sideExtremeValue(const std::vector<int> & boundary_id, const field::NekFieldEnum & field,
614 const nek_mesh::NekMeshEnum pp_mesh, const bool max);
615
620int Nfaces();
621
627bool isHeatFluxBoundary(const int boundary);
628
634bool isMovingMeshBoundary(const int boundary);
635
641bool isTemperatureBoundary(const int boundary);
642
648const std::string temperatureBoundaryType(const int boundary);
649
655
661
666int dim();
667
676
682
688
696bool
697validBoundaryIDs(const std::vector<int> & boundary_id, int & first_invalid_id, int & n_boundaries);
698
704void storeBoundaryCoupling(const std::vector<int> & boundary_id, int & N);
705
711{
713 int flux = -1;
714
716 int heat_source = -1;
717
720
723
726
729
732
735
738
741};
742
748{
749 double U_ref = 1;
750 double T_ref = 0;
751 double dT_ref = 1;
752 double P_ref = 1;
753 double L_ref = 1;
754 double A_ref = 1;
755 double V_ref = 1;
756 double rho_ref = 1;
757 double Cp_ref = 1;
758 double flux_ref = 1;
759 double source_ref = 1;
760 double t_ref = 1;
761 double s01_ref = 0;
762 double ds01_ref = 1;
763 double s02_ref = 0;
764 double ds02_ref = 1;
765 double s03_ref = 0;
766 double ds03_ref = 1;
767};
768
775
780void (*solutionPointer(const field::NekWriteEnum & field))(int, dfloat);
781
787double scalar01(const int id);
788
794double scalar02(const int id);
795
801double scalar03(const int id);
802
808double usrwrk00(const int id);
809
815double usrwrk01(const int id);
816
822double usrwrk02(const int id);
823
829double temperature(const int id);
830
836double pressure(const int id);
837
843double unity(const int id);
844
850double velocity_x(const int id);
851
857double velocity_y(const int id);
858
864double velocity_z(const int id);
865
871double velocity(const int id);
872
878double velocity_x_squared(const int id);
879
885double velocity_y_squared(const int id);
886
892double velocity_z_squared(const int id);
893
899void flux(const int id, const dfloat value);
900
906void heat_source(const int id, const dfloat value);
907
913void x_displacement(const int id, const dfloat value);
914
920void y_displacement(const int id, const dfloat value);
921
927void z_displacement(const int id, const dfloat value);
928
938void initializeDimensionalScales(const double U,
939 const double T,
940 const double dT,
941 const double L,
942 const double rho,
943 const double Cp,
944 const double s01,
945 const double ds01,
946 const double s02,
947 const double ds02,
948 const double s03,
949 const double ds03);
950
964void dimensionalize(const field::NekFieldEnum & field, double & value);
965
971
977
983
989
995
1001
1007
1013
1014// useful concept from Stack Overflow for templating MPI calls
1015template <typename T>
1016MPI_Datatype resolveType();
1017
1025template <typename T>
1026void
1027allgatherv(const std::vector<int> & base_counts, const T * input, T * output, const int multiplier)
1028{
1029 int * recvCounts = (int *)calloc(commSize(), sizeof(int));
1030 int * displacement = (int *)calloc(commSize(), sizeof(int));
1031 displacementAndCounts(base_counts, recvCounts, displacement, multiplier);
1032
1033 MPI_Allgatherv(input,
1034 recvCounts[commRank()],
1036 output,
1037 (const int *)recvCounts,
1038 (const int *)displacement,
1040 platform->comm.mpiComm);
1041
1044}
1045
1046} // end namespace nekrs
Definition NekBoundaryCoupling.h:30
Definition CardinalEnums.h:89
NekFieldEnum
Enumeration of possible fields to read from nekRS.
Definition CardinalEnums.h:92
NekWriteEnum
Enumeration of possible fields to write in nekRS.
Definition CardinalEnums.h:114
NekMeshEnum
Definition CardinalEnums.h:48
Cardinal-specific nekRS API.
Definition NekUility.C:22
int NboundaryID()
double usrwrk02(const int id)
void write_field_file(const std::string &prefix, const dfloat time, const int &step)
void limitTemperature(const double *min_T, const double *max_T)
double volumeExtremeValue(const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
void interpolateVolumeHex3D(const double *I, double *x, int N, double *Ix, int M)
double referenceTime()
Real scratchUnits(const int slot)
bool hasTemperatureVariable()
double referenceVelocity()
double velocity_x(const int id)
void allgatherv(const std::vector< int > &base_counts, const T *input, T *output, const int multiplier=1)
Definition NekInterface.h:1027
Point gllPoint(int local_elem_id, int local_node_id)
void dimensionalizeSideIntegral(const field::NekFieldEnum &integrand, const Real &area, double &integral)
dfloat * getVgeo()
double velocity_y(const int id)
int polynomialOrder()
Real referenceAdditiveScale(const field::NekFieldEnum &field)
int buildOnly()
void setAbsoluteTol(double tol)
void storeBoundaryCoupling(const std::vector< int > &boundary_id, int &N)
void setStartTime(const double &start)
double scalar03(const int id)
void write_usrwrk_field_file(const int &slot, const std::string &prefix, const dfloat &time, const int &step, const bool &write_coords)
Point gllPointFace(int local_elem_id, int local_face_id, int local_node_id)
void gradient(const int offset, const double *f, double *grad_f, const nek_mesh::NekMeshEnum pp_mesh)
void dimensionalize(const field::NekFieldEnum &field, double &value)
Dimensionalize a field by multiplying the nondimensional form by the reference.
void interpolateSurfaceFaceHex3D(double *scratch, const double *I, double *x, int N, double *Ix, int M)
double referencePressure()
bool hasTemperatureSolve()
double referenceFlux()
bool endControlNumSteps()
int velocityFieldOffset()
mesh_t * flowMesh()
mesh_t * getMesh(const nek_mesh::NekMeshEnum pp_mesh)
void initializeDimensionalScales(const double U, const double T, const double dT, const double L, const double rho, const double Cp, const double s01, const double ds01, const double s02, const double ds02, const double s03, const double ds03)
void setNekSetupTime(const double &time)
double referenceArea()
void heat_source(const int id, const dfloat value)
int Nelements()
double referenceLength()
double sideMassFluxWeightedIntegral(const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
bool endControlElapsedTime()
double referenceVolume()
int commRank()
double Pr()
double sideIntegral(const std::vector< int > &boundary_id, const field::NekFieldEnum &integrand, const nek_mesh::NekMeshEnum pp_mesh)
int fieldOffset()
int Nfaces()
double velocity_x_squared(const int id)
void dimensionalizeVolumeIntegral(const field::NekFieldEnum &integrand, const Real &volume, double &integral)
double volumeIntegral(const field::NekFieldEnum &integrand, const double &volume, const nek_mesh::NekMeshEnum pp_mesh)
double viscosity()
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)
double scalar02(const int id)
void dimensionalizeArea(double &integral)
const std::string temperatureBoundaryType(const int boundary)
void scaleUsrwrk(const unsigned int &slot, const dfloat &value)
int NfaceVertices()
Number of vertices required to define an element face Vertices refer to the points required to place ...
void x_displacement(const int id, const dfloat value)
bool hasUserMeshSolver()
void interpolationMatrix(double *I, int starting_points, int ending_points)
double(*)(int) solutionPointer(const field::NekFieldEnum &field)
Definition NekInterface.h:774
bool isHeatFluxBoundary(const int boundary)
int commSize()
double unity(const int id)
double temperature(const int id)
MPI_Datatype resolveType()
std::vector< double > usrwrkSideIntegral(const unsigned int &slot, const std::vector< int > &boundary, const nek_mesh::NekMeshEnum pp_mesh)
void initializeScratch(const unsigned int &n_slots)
void updateHostMeshParameters()
Update the mesh parameters on host.
double heatFluxIntegral(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
double area(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
Point centroidFace(int local_elem_id, int local_face_id)
bool normalizeFlux(const NekBoundaryCoupling &nek_boundary_coupling, const std::vector< int > &boundary, const double moose_integral, double nek_integral, double &normalized_nek_integral)
Point centroid(int local_elem_id)
bool hasBlendingSolver()
int NboundaryFaces()
bool hasCHT()
mesh_t * entireMesh()
void dimensionalizeVolume(double &integral)
bool endControlTime()
dfloat * getSgeo()
bool isInitialized()
bool hasVariableDt()
bool isTemperatureBoundary(const int boundary)
double getNekSetupTime()
double pressure(const int id)
bool validBoundaryIDs(const std::vector< int > &boundary_id, int &first_invalid_id, int &n_boundaries)
bool hasScalarVariable(int scalarId)
void copyDeformationToDevice()
Copy the deformation from host to device.
double massFlowrate(const std::vector< int > &boundary_id, const nek_mesh::NekMeshEnum pp_mesh)
void checkFieldValidity(const field::NekFieldEnum &field)
double usrwrkVolumeIntegral(const unsigned int &slot, const nek_mesh::NekMeshEnum pp_mesh)
bool hasHeatSourceKernel()
double volume(const nek_mesh::NekMeshEnum pp_mesh)
double usrwrk00(const int id)
void nondimensional(const bool n)
int scalarFieldOffset()
bool isMovingMeshBoundary(const int boundary)
double velocity_z(const int id)
void setRelativeTol(double tol)
double referenceSource()
double usrwrk01(const int id)
void freeScratch()
Free the scratch space.
double pressureSurfaceForce(const std::vector< int > &boundary_id, const Point &direction, const nek_mesh::NekMeshEnum pp_mesh)
void y_displacement(const int id, const dfloat value)
double sideExtremeValue(const std::vector< int > &boundary_id, const field::NekFieldEnum &field, const nek_mesh::NekMeshEnum pp_mesh, const bool max)
int dim()
double velocity_z_squared(const int id)
void initializeHostMeshParameters()
Allocate memory for the host mesh parameters.
double velocity(const int id)
bool hasMovingMesh()
bool scratchAvailable()
mesh_t * temperatureMesh()
void flux(const int id, const dfloat value)
double scalar01(const int id)
void z_displacement(const int id, const dfloat value)
void displacementAndCounts(const std::vector< int > &base_counts, int *counts, int *displacement, const int multiplier)
double velocity_y_squared(const int id)
Definition NekInterface.h:748
double A_ref
Definition NekInterface.h:754
double flux_ref
Definition NekInterface.h:758
double s02_ref
Definition NekInterface.h:763
double source_ref
Definition NekInterface.h:759
double Cp_ref
Definition NekInterface.h:757
double dT_ref
Definition NekInterface.h:751
double ds03_ref
Definition NekInterface.h:766
double T_ref
Definition NekInterface.h:750
double P_ref
Definition NekInterface.h:752
double V_ref
Definition NekInterface.h:755
double L_ref
Definition NekInterface.h:753
double t_ref
Definition NekInterface.h:760
double s01_ref
Definition NekInterface.h:761
double ds02_ref
Definition NekInterface.h:764
double ds01_ref
Definition NekInterface.h:762
double rho_ref
Definition NekInterface.h:756
double s03_ref
Definition NekInterface.h:765
double U_ref
Definition NekInterface.h:749
Definition NekInterface.h:711
int boundary_temperature
boundary temperature (for separate domain coupling)
Definition NekInterface.h:731
int mesh_velocity_z
z-velocity of moving boundary (for mesh blending solver)
Definition NekInterface.h:725
int flux
boundary heat flux (for conjugate heat transfer)
Definition NekInterface.h:713
int mesh_velocity_y
y-velocity of moving boundary (for mesh blending solver)
Definition NekInterface.h:722
int heat_source
volumetric heat source (for volumetric heating)
Definition NekInterface.h:716
int boundary_velocity
boundary velocity (for separate domain coupling)
Definition NekInterface.h:728
int boundary_scalar03
boundary scalar03 (for separate domain coupling)
Definition NekInterface.h:740
int boundary_scalar01
boundary scalar01 (for separate domain coupling)
Definition NekInterface.h:734
int mesh_velocity_x
x-velocity of moving boundary (for mesh blending solver)
Definition NekInterface.h:719
int boundary_scalar02
boundary scalar02 (for separate domain coupling)
Definition NekInterface.h:737