21 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_
22 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_
26 #include "finiteElement/PDEUtilities.hpp"
27 #include "LvArray/src/tensorOps.hpp"
38 template<
int NUM_SUPPORT_POINTS,
40 int NUM_QUADRATURE_POINTS >
68 #pragma nv_diag_suppress 20012
106 template<
typename SUBREGION_TYPE >
125 template<
typename STACK_VARIABLES_TYPE >
148 template<
typename STACK_VARIABLES_TYPE >
175 real64 (& samplingPointCoord)[3] )
191 template<
typename SUBREGION_TYPE,
192 typename MESH_DATA_TYPE >
196 SUBREGION_TYPE
const & cellSubRegion,
197 MESH_DATA_TYPE & meshData )
216 template<
typename LEAF,
217 typename SUBREGION_TYPE,
218 typename MESH_DATA_TYPE >
222 SUBREGION_TYPE
const & cellSubRegion,
223 MESH_DATA_TYPE & meshData
226 LEAF::template fillMeshData< SUBREGION_TYPE >( nodeManager,
240 template<
typename MESH_DATA_TYPE,
241 typename STACK_VARIABLES_TYPE >
245 MESH_DATA_TYPE
const & meshData,
246 STACK_VARIABLES_TYPE & stack )
263 template<
typename LEAF,
264 typename MESH_DATA_TYPE >
267 MESH_DATA_TYPE
const & meshData,
268 typename LEAF::StackVariables & stack )
const
270 LEAF::setupStack( cellIndex, meshData, stack );
284 if constexpr ( N == 1 )
286 return PDEUtilities::FunctionSpace::H1;
288 else if constexpr ( N == 3 )
290 return PDEUtilities::FunctionSpace::H1vector;
294 static_assert( N == 1 || N == 3,
"Unsupported number of components per support point" );
310 template<
localIndex NUMDOFSPERTRIALSUPPORTPOINT,
313 typename STACK_VARIABLES_TYPE >
317 real64 ( & matrix )[MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
318 real64 const & scaleFactor )
336 template<
typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT,
bool UPPER = false >
339 real64 ( & matrix )[LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
340 real64 const scaleFactor = 1.0 )
const
343 LEAF::maxSupportPoints,
362 template<
localIndex NUMDOFSPERTRIALSUPPORTPOINT,
364 typename STACK_VARIABLES_TYPE >
368 real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
369 real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
370 real64 const scaleFactor )
391 template<
typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT >
396 real64 const ( &dofs )[LEAF::maxSupportPoints]
397 [NUMDOFSPERTRIALSUPPORTPOINT],
398 real64 ( & targetVector )[LEAF::maxSupportPoints]
399 [NUMDOFSPERTRIALSUPPORTPOINT],
400 real64 const scaleFactor = 1.0 )
const
403 addEvaluatedGradGradStabilization< NUMDOFSPERTRIALSUPPORTPOINT, LEAF::maxSupportPoints >( stack,
429 m_numSupportPoints( numSupportPoints ),
430 m_maxSupportPoints( maxSupportPoints ),
431 m_numQuadraturePoints( numQuadraturePoints )
445 return m_numQuadraturePoints;
#define GEOS_HOST_DEVICE
Marks a host-device function.
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
#define GEOS_ERROR(...)
Raise a hard error and terminate the program.
This class provides an interface to ObjectManagerBase in order to manage edge data.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Device-compatible (non virtual) Base class for all finite element formulations.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(STACK_VARIABLES_TYPE const &stack)
Get the number of quadrature points.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints()
Get the number of quadrature points.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addGradGradStabilization(STACK_VARIABLES_TYPE const &stack, real64(&matrix)[MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT], real64 const &scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
GEOS_HOST_DEVICE ~FiniteElementBase_impl()=default
Default destructor.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
GEOS_HOST_DEVICE FiniteElementBase_impl(FiniteElementBase_impl const &)=default
Default copy constructor.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(STACK_VARIABLES_TYPE const &stack)
Get the number of support points.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const &linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
static void fillMeshData(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MESH_DATA_TYPE &meshData)
Method to fill a MeshData object.
GEOS_HOST_DEVICE void setup(localIndex const &cellIndex, MESH_DATA_TYPE const &meshData, typename LEAF::StackVariables &stack) const
Abstract setup method, possibly computing cell-dependent properties.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addEvaluatedGradGradStabilization(STACK_VARIABLES_TYPE const &stack, real64 const (&dofs)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
GEOS_HOST_DEVICE FiniteElementBase_impl()=default
Default constructor.
constexpr static localIndex numSupportPoints
The number of support points per element.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl &&)=default
Default move assignment operator.
static GEOS_HOST_DEVICE localIndex getMaxSupportPoints()
Get the maximum number of support points.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
GEOS_HOST_DEVICE void addGradGradStabilizationMatrix(typename LEAF::StackVariables const &stack, real64(&matrix)[LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add stabilization of grad-grad bilinear form to input matrix.
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void addEvaluatedGradGradStabilizationVector(typename LEAF::StackVariables const &stack, real64 const (&dofs)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add a grad-grad stabilization operator evaluated at a provided vector of dofs to input vector.
constexpr static localIndex numNodes
The number of nodes per element.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints()
Get the number of support points.
constexpr static GEOS_HOST_DEVICE PDEUtilities::FunctionSpace getFunctionSpace()
Getter for the function space.
static void initialize(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MESH_DATA_TYPE &meshData)
Abstract initialization method.
constexpr static int numSamplingPoints
The number of sampling points per element.
GEOS_HOST_DEVICE FiniteElementBase_impl(FiniteElementBase_impl &&)=default
Default move constructor.
constexpr static localIndex numFaces
The number of faces per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void setupStack(localIndex const &cellIndex, MESH_DATA_TYPE const &meshData, STACK_VARIABLES_TYPE &stack)
Empty setup method.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl const &)=default
Default copy assignment operator.
Base class for FEM element implementations.
localIndex getNumSupportPoints() const
Getter for the number of support points per element.
virtual ~FiniteElementBase()=default
Destructor.
localIndex getNumQuadraturePoints() const
Getter for the number of quadrature points per element.
localIndex getMaxSupportPoints() const
Get the maximum number of support points for this element.
FiniteElementBase(localIndex const numSupportPoints, localIndex const maxSupportPoints, localIndex const numQuadraturePoints)
Default constructor.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Variables used to initialize the class.
Kernel variables allocated on the stack.