20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TRIANGLEFACELAGRANGE1GAUSS_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TRIANGLEFACELAGRANGE1GAUSS_HPP_
50 template<
typename NUM_Q_POINTS >
56 static_assert( ( NUM_Q_POINTS::value == 1 ||
57 NUM_Q_POINTS::value == 4 ||
58 NUM_Q_POINTS::value == 6 ),
59 "NUM_Q_POINTS::value must be 1, 4, or 6!" );
165 real64 const r = pointCoord[0];
166 real64 const s = pointCoord[1];
168 N[0] = (1.0 - r - s) * r * s;
184 real64 const qCoords[2] = {quadratureParentCoords0( q ), quadratureParentCoords1( q ) };
209 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT,
bool UPPER >
216 real64 const & scaleFactor );
234 constexpr
static real64 parentArea = 0.5;
238 constexpr
static real64 weight = parentArea;
355 template<
typename NUM_Q_POINTS >
356 template< localIndex NUMDOFSPERTRIALSUPPORTPOINT,
bool UPPER >
362 [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
363 [maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
364 real64 const & scaleFactor )
371 template<
typename NUM_Q_POINTS >
377 real64 ( & N )[numNodes] )
379 real64 const r = coords[0];
380 real64 const s = coords[1];
387 template<
typename NUM_Q_POINTS >
396 real64 const qCoords[2] = {quadratureParentCoords0( q ), quadratureParentCoords1( q ) };
402 template<
typename NUM_Q_POINTS >
408 real64 ( & N )[numNodes] )
410 return calcN( q, N );
415 template<
typename NUM_Q_POINTS >
421 real64 const (&X)[numNodes][3] )
424 real64 n[3] = { ( X[1][1] - X[0][1] ) * ( X[2][2] - X[0][2] ) - ( X[2][1] - X[0][1] ) * ( X[1][2] - X[0][2] ),
425 ( X[2][0] - X[0][0] ) * ( X[1][2] - X[0][2] ) - ( X[1][0] - X[0][0] ) * ( X[2][2] - X[0][2] ),
426 ( X[1][0] - X[0][0] ) * ( X[2][1] - X[0][1] ) - ( X[2][0] - X[0][0] ) * ( X[1][1] - X[0][1] )};
428 return sqrt( n[0] * n[0] + n[1] * n[1] + n[2] * n[2] ) * weight * quadratureWeight( q );
#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_UNUSED_PARAM(X)
Mark an unused argument and silence compiler warnings.
Base class for FEM element implementations.
static GEOS_HOST_DEVICE void calcN(localIndex const q, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature point.
static GEOS_HOST_DEVICE void calcBubbleN(localIndex const q, real64(&N)[1])
Calculate shape bubble functions values at a quadrature point.
constexpr static localIndex numNodes
The number of nodes/support points per element.
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const override
Virtual getter for the number of quadrature points per element.
static GEOS_HOST_DEVICE void calcN(localIndex const q, StackVariables const &stack, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature point.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature point.
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const override
Get the maximum number of support points for this element.
static GEOS_HOST_DEVICE void getPermutation(int(&permutation)[numNodes])
Calculate the node permutation between the parent element and the geometric element....
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN(real64 const (&coords)[2], real64(&N)[numNodes])
Calculate shape functions values at a single point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcBubbleN(real64 const (&pointCoord)[2], real64(&N)[1])
Calculate shape bubble functions values at a given point in the parent space.
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
static GEOS_HOST_DEVICE void addGradGradStabilization(StackVariables 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...
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Kernel variables allocated on the stack.