20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_HPP_
132 real64 (& samplingPointCoord)[3] )
140 real64 const r = i0 * step;
141 real64 const s = i1 * step;
142 real64 const t = i2 * 2 * step;
145 samplingPointCoord[0] = r;
146 samplingPointCoord[1] = s;
147 samplingPointCoord[2] = t;
152 samplingPointCoord[0] = 1 - r;
153 samplingPointCoord[1] = 1 - s;
154 samplingPointCoord[2] = t;
207 GEOS_ERROR(
"Unsupported bubble functions for wedge elements" );
223 GEOS_ERROR(
"Unsupported bubble functions for wedge elements" );
309 jacobianTransformation( q, X, J );
310 return LvArray::tensorOps::invert< 3 >( J );
316 constexpr
static real64 parentVolume = 1.0;
323 constexpr
static real64 quadratureCrossSectionCoord = 1.0 / 6.0;
328 constexpr
static real64 quadratureLongitudinalCoord = 1.0 / 1.732050807568877293528;
337 template<
typename T >
340 constexpr
static T linearMap( T
const indexT, T
const indexL )
342 return 2 * indexT + indexL;
355 return 0.5 * ( a & 2 );
368 return 0.25 * ( a & 4 );
381 return -1.0 + 2 * ( a & 1 );
394 return quadratureCrossSectionCoord + 0.5 * parentCoords0( q );
407 return quadratureCrossSectionCoord + 0.5 * parentCoords1( q );
420 return parentCoords2( q ) * quadratureLongitudinalCoord;
430 static void jacobianTransformation(
int const q,
445 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
446 real64 const ( &invJ )[3][3],
456 H1_Wedge_Lagrange1_Gauss6::
457 jacobianTransformation(
int const q,
458 real64 const (&X)[numNodes][3],
461 real64 const r = quadratureParentCoords0( q );
462 real64 const s = quadratureParentCoords1( q );
463 real64 const xi = quadratureParentCoords2( q );
465 real64 const psiTRI[3] = { 1.0 - r - s, r, s };
466 real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
467 constexpr
real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
468 constexpr
real64 dpsiLIN[2] = { -0.5, 0.5 };
474 real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
475 dpsiTRI[1][a] * psiLIN[b],
476 psiTRI[a] * dpsiLIN[b] };
477 localIndex const nodeIndex = linearMap( a, b );
478 for(
int i = 0; i < 3; ++i )
480 for(
int j = 0; j < 3; ++j )
482 J[i][j] = J[i][j] + dNdXi[ j ] * X[nodeIndex][i];
494 H1_Wedge_Lagrange1_Gauss6::
495 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
496 real64 const ( &invJ )[3][3],
497 real64 (& gradN)[numNodes][3] )
499 real64 const r = quadratureParentCoords0( q );
500 real64 const s = quadratureParentCoords1( q );
501 real64 const xi = quadratureParentCoords2( q );
503 real64 const psiTRI[3] = { 1.0 - r - s, r, s };
504 real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
505 constexpr
real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
506 constexpr
real64 dpsiLIN[2] = { -0.5, 0.5 };
512 real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
513 dpsiTRI[1][a] * psiLIN[b],
514 psiTRI[a] * dpsiLIN[b] };
515 localIndex const nodeIndex = linearMap( a, b );
516 for(
int i = 0; i < 3; ++i )
518 gradN[nodeIndex][i] = 0.0;
519 for(
int j = 0; j < 3; ++j )
521 gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
536 real64 const r = coords[0];
537 real64 const s = coords[1];
538 real64 const xi = coords[2];
540 N[0] = 0.5*( 1.0 - r - s ) * ( 1.0 - xi );
541 N[1] = 0.5*( 1.0 - r - s ) * ( 1.0 + xi );
542 N[2] = 0.5* r * ( 1.0 - xi );
543 N[3] = 0.5* r * ( 1.0 + xi );
544 N[4] = 0.5* s * ( 1.0 - xi );
545 N[5] = 0.5* s * ( 1.0 + xi );
557 real64 const pointCoord[3] = {quadratureParentCoords0( q ),
558 quadratureParentCoords1( q ),
559 quadratureParentCoords2( q )};
561 calcN( pointCoord, N );
569 real64 ( & N )[numNodes] )
571 return calcN( q, N );
581 real64 const (&X)[numNodes][3],
582 real64 (& gradN)[numNodes][3] )
586 jacobianTransformation( q, X, J );
588 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
590 applyJacobianTransformationToShapeFunctionsDerivatives( q, J, gradN );
592 return detJ * weight;
599 real64 const (&X)[numNodes][3],
601 real64 ( & gradN )[numNodes][3] )
610 real64 const (&X)[numNodes][3],
611 real64 (& gradN)[numFaces][3] )
614 GEOS_ERROR(
"Unsupported bubble functions for wedge elements" );
625 real64 const (&X)[numNodes][3] )
629 jacobianTransformation( q, X, J );
631 return LvArray::tensorOps::determinant< 3 >( J ) * weight;
#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.
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Base class for FEM element implementations.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate shape functions values for each support face at a quadrature point.
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const &linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
constexpr static localIndex numNodes
The number of nodes/support points per element.
constexpr static localIndex numFaces
The number of faces/support points per element.
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 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the shape bubble function derivatives wrt the physical coordinates.
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.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
constexpr static int numSamplingPoints
The number of sampling 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.
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 calcN(real64 const (&coords)[3], real64(&N)[numNodes])
Calculate shape functions values at a single point.
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack)
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.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
static GEOS_HOST_DEVICE real64 invJacobianTransformation(int const q, real64 const (&X)[numNodes][3], real64(&J)[3][3])
Calculates the isoparametric "Jacobian" transformation matrix/mapping from the parent space to the ph...
static GEOS_HOST_DEVICE real64 calcGradN(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack, real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcFaceBubbleN(real64 const (&pointCoord)[3], real64(&N)[numFaces])
Calculate shape functions values for each support face at a given point in the parent space.
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.