20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1WEDGELAGRANGE1GAUSS6_HPP_
69 template<
typename SubregionType >
114 real64 (& samplingPointCoord)[3] )
122 real64 const r = i0 * step;
123 real64 const s = i1 * step;
124 real64 const t = i2 * 2 * step;
127 samplingPointCoord[0] = r;
128 samplingPointCoord[1] = s;
129 samplingPointCoord[2] = t;
134 samplingPointCoord[0] = 1 - r;
135 samplingPointCoord[1] = 1 - s;
136 samplingPointCoord[2] = t;
189 real64 const r = pointCoord[0];
190 real64 const s = pointCoord[1];
191 real64 const xi = pointCoord[2];
214 real64 const pointCoord[3] = {quadratureParentCoords0( q ),
215 quadratureParentCoords1( q ),
216 quadratureParentCoords2( q )};
334 jacobianTransformation( q, X, J );
335 return LvArray::tensorOps::invert< 3 >( J );
341 constexpr
static real64 parentVolume = 1.0;
348 constexpr
static real64 quadratureCrossSectionCoord = 1.0 / 6.0;
353 constexpr
static real64 quadratureLongitudinalCoord = 1.0 / 1.732050807568877293528;
362 template<
typename T >
365 constexpr
static T linearMap( T
const indexT, T
const indexL )
367 return 2 * indexT + indexL;
380 return 0.5 * ( a & 2 );
393 return 0.25 * ( a & 4 );
406 return -1.0 + 2 * ( a & 1 );
419 return quadratureCrossSectionCoord + 0.5 * parentCoords0( q );
432 return quadratureCrossSectionCoord + 0.5 * parentCoords1( q );
445 return parentCoords2( q ) * quadratureLongitudinalCoord;
455 static void jacobianTransformation(
int const q,
470 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
471 real64 const ( &invJ )[3][3],
481 H1_Wedge_Lagrange1_Gauss6_impl::
482 jacobianTransformation(
int const q,
483 real64 const (&X)[numNodes][3],
486 real64 const r = quadratureParentCoords0( q );
487 real64 const s = quadratureParentCoords1( q );
488 real64 const xi = quadratureParentCoords2( q );
490 real64 const psiTRI[3] = { 1.0 - r - s, r, s };
491 real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
492 constexpr
real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
493 constexpr
real64 dpsiLIN[2] = { -0.5, 0.5 };
499 real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
500 dpsiTRI[1][a] * psiLIN[b],
501 psiTRI[a] * dpsiLIN[b] };
502 localIndex const nodeIndex = linearMap( a, b );
503 for(
int i = 0; i < 3; ++i )
505 for(
int j = 0; j < 3; ++j )
507 J[i][j] = J[i][j] + dNdXi[ j ] * X[nodeIndex][i];
519 H1_Wedge_Lagrange1_Gauss6_impl::
520 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
521 real64 const ( &invJ )[3][3],
522 real64 (& gradN)[numNodes][3] )
524 real64 const r = quadratureParentCoords0( q );
525 real64 const s = quadratureParentCoords1( q );
526 real64 const xi = quadratureParentCoords2( q );
528 real64 const psiTRI[3] = { 1.0 - r - s, r, s };
529 real64 const psiLIN[2] = { 0.5 - 0.5*xi, 0.5 + 0.5*xi };
530 constexpr
real64 dpsiTRI[2][3] = { { -1.0, 1.0, 0.0 }, { -1.0, 0.0, 1.0 } };
531 constexpr
real64 dpsiLIN[2] = { -0.5, 0.5 };
537 real64 const dNdXi[3] = { dpsiTRI[0][a] * psiLIN[b],
538 dpsiTRI[1][a] * psiLIN[b],
539 psiTRI[a] * dpsiLIN[b] };
540 localIndex const nodeIndex = linearMap( a, b );
541 for(
int i = 0; i < 3; ++i )
543 gradN[nodeIndex][i] = 0.0;
544 for(
int j = 0; j < 3; ++j )
546 gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
561 real64 const r = coords[0];
562 real64 const s = coords[1];
563 real64 const xi = coords[2];
565 N[0] = 0.5*( 1.0 - r - s ) * ( 1.0 - xi );
566 N[1] = 0.5*( 1.0 - r - s ) * ( 1.0 + xi );
567 N[2] = 0.5* r * ( 1.0 - xi );
568 N[3] = 0.5* r * ( 1.0 + xi );
569 N[4] = 0.5* s * ( 1.0 - xi );
570 N[5] = 0.5* s * ( 1.0 + xi );
582 real64 const pointCoord[3] = {quadratureParentCoords0( q ),
583 quadratureParentCoords1( q ),
584 quadratureParentCoords2( q )};
586 calcN( pointCoord, N );
594 real64 ( & N )[numNodes] )
596 return calcN( q, N );
605 real64 const (&X)[numNodes][3],
608 for(
int i = 0; i < 3; ++i )
610 for(
int j = 0; j < 3; ++j )
615 jacobianTransformation( q, X, J );
616 real64 const detJ = LvArray::tensorOps::determinant< 3 >( J );
617 return detJ * weight;
624 real64 const (&X)[numNodes][3],
636 real64 const (&X)[numNodes][3],
637 real64 (& gradN)[numNodes][3] )
641 jacobianTransformation( q, X, J );
643 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
645 applyJacobianTransformationToShapeFunctionsDerivatives( q, J, gradN );
647 return detJ * weight;
654 real64 const (&X)[numNodes][3],
656 real64 ( & gradN )[numNodes][3] )
665 real64 const (&X)[numNodes][3],
666 real64 (& gradN)[numFaces][3] )
671 jacobianTransformation( q, X, J );
673 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
677 real64 const r = quadratureParentCoords0( q );
678 real64 const s = quadratureParentCoords1( q );
679 real64 const xi = quadratureParentCoords2( q );
703 gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
704 gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
705 gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
708 return detJ * weight;
718 real64 const (&X)[numNodes][3] )
722 jacobianTransformation( q, X, J );
724 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.
Device-compatible (non virtual) Base class for all finite element formulations.
constexpr static localIndex numQuadraturePoints
The number of quadrature points per element.
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
constexpr static localIndex maxSupportPoints
The maximum number of support points per element.
constexpr static localIndex numNodes
The number of nodes per element.
constexpr static localIndex numFaces
The number of faces per element.
GEOS_HOST_DEVICE FiniteElementBase_impl & operator=(FiniteElementBase_impl const &)=default
Default copy assignment operator.
Base class for FEM element implementations.
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.
DO_NOT_DOCUMENT static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
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 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 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the shape bubble function derivatives wrt the physical coordinates.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 calcJacobian(localIndex const q, real64 const (&X)[numNodes][3], real64(&J)[3][3])
Calculate the Jacobian matrix at a quadrature point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const &linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 calcJacobian(localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack, real64(&J)[3][3])
Calculate the Jacobian matrix at a quadrature point.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate face bubble functions values for each face at a quadrature point.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
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.
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.
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 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 real64 calcGradN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numNodes][3])
Calculate the shape functions derivatives wrt the physical coordinates.
H1_Wedge_Lagrange1_Gauss6_impl * getImpl()
Get the device-compatible implementation type.
H1_Wedge_Lagrange1_Gauss6_impl const * getImpl() const
Get the device-compatible implementation type.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 gradientBubble(const real64 xi)
The gradient of the bubble basis function for support point 1 evaluated at a point along the axes.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 valueBubble(const real64 xi)
The value of the bubble basis function.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 gradient(const int index, const real64 xi)
The gradient of the basis function for a support point evaluated at a point along the axes.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 value(const int index, const real64 xi)
The value of the basis function for a support point evaluated at a point along the axes.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Kernel variables (dof numbers, jacobian and residual) located on the stack.