20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TETRAHEDRONLAGRANGE1GAUSS_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1TETRAHEDRONLAGRANGE1GAUSS_HPP_
49 template<
typename NUM_Q_POINTS >
60 template<
typename SubregionType >
85 static_assert( ( NUM_Q_POINTS::value == 1 ||
86 NUM_Q_POINTS::value == 5 ||
87 NUM_Q_POINTS::value == 14 ),
88 "NUM_Q_POINTS::value must be 1, 5, or 14!" );
157 real64 (& samplingPointCoord)[3] )
165 real64 const r = i0 * step;
166 real64 const s = i1 * step;
167 real64 const t = i2 * step;
170 samplingPointCoord[0] = r;
171 samplingPointCoord[1] = s;
172 samplingPointCoord[2] = t;
177 samplingPointCoord[0] = 1 - t - r;
178 samplingPointCoord[1] = 1 - t - s;
179 samplingPointCoord[2] = t;
233 real64 const r = pointCoord[0];
234 real64 const s = pointCoord[1];
235 real64 const t = pointCoord[2];
237 N[0] = (1 - r - s - t) * r * t;
238 N[1] = (1 - r - s - t) * r * s;
239 N[2] = (1 - r - s - t) * t * s;
256 real64 const pointCoord[3] = {quadratureParentCoords0( q ),
257 quadratureParentCoords1( q ),
258 quadratureParentCoords2( q )};
378 return LvArray::tensorOps::invert< 3 >( J );
382 constexpr
static real64 parentVolume = 1.0 / 6.0;
385 constexpr
static real64 weight = parentVolume;
410 0.073493043116361949544,
411 0.073493043116361949544,
412 0.073493043116361949544,
413 0.11268792571801585080,
414 0.11268792571801585080,
415 0.11268792571801585080,
416 0.11268792571801585080,
417 0.042546020777081466438,
418 0.042546020777081466438,
419 0.042546020777081466438,
420 0.042546020777081466438,
421 0.042546020777081466438,
422 0.042546020777081466438 };
452 0.092735250310891226402,
453 0.092735250310891226402,
454 0.092735250310891226402,
455 0.067342242210098170608,
456 0.31088591926330060980,
457 0.31088591926330060980,
458 0.31088591926330060980,
459 0.045503704125649649492,
460 0.045503704125649649492,
461 0.045503704125649649492,
462 0.45449629587435035051,
463 0.45449629587435035051,
464 0.45449629587435035051 };
494 0.72179424906732632079,
495 0.092735250310891226402,
496 0.092735250310891226402,
497 0.31088591926330060980,
498 0.067342242210098170608,
499 0.31088591926330060980,
500 0.31088591926330060980,
501 0.045503704125649649492,
502 0.45449629587435035051,
503 0.45449629587435035051,
504 0.045503704125649649492,
505 0.045503704125649649492,
506 0.45449629587435035051 };
536 0.092735250310891226402,
537 0.72179424906732632079,
538 0.092735250310891226402,
539 0.31088591926330060980,
540 0.31088591926330060980,
541 0.067342242210098170608,
542 0.31088591926330060980,
543 0.45449629587435035051,
544 0.045503704125649649492,
545 0.45449629587435035051,
546 0.045503704125649649492,
547 0.45449629587435035051,
548 0.045503704125649649492 };
568 template<
typename NUM_Q_POINTS >
572 H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS >::
573 determinantJacobianTransformation(
real64 const (&X)[numNodes][3] )
575 return ( X[1][0] - X[0][0] )*( ( X[2][1] - X[0][1] )*( X[3][2] - X[0][2] ) - ( X[3][1] - X[0][1] )*( X[2][2] - X[0][2] ) )
576 + ( X[1][1] - X[0][1] )*( ( X[3][0] - X[0][0] )*( X[2][2] - X[0][2] ) - ( X[2][0] - X[0][0] )*( X[3][2] - X[0][2] ) )
577 + ( X[1][2] - X[0][2] )*( ( X[2][0] - X[0][0] )*( X[3][1] - X[0][1] ) - ( X[3][0] - X[0][0] )*( X[2][1] - X[0][1] ) );
582 template<
typename NUM_Q_POINTS >
590 real64 const r = coords[0];
591 real64 const s = coords[1];
592 real64 const t = coords[2];
595 N[0] = 1 - r - s - t;
602 template<
typename NUM_Q_POINTS >
610 real64 const pointCoord[3] = {quadratureParentCoords0( q ),
611 quadratureParentCoords1( q ),
612 quadratureParentCoords2( q )};
614 calcN( pointCoord, N );
617 template<
typename NUM_Q_POINTS >
623 real64 ( & N )[numNodes] )
625 return calcN( q, N );
629 template<
typename NUM_Q_POINTS >
634 real64 const (&X)[numNodes][3],
637 for(
int i = 0; i < 3; ++i )
639 for(
int j = 0; j < 3; ++j )
644 J[0][0] = X[1][0] - X[0][0];
645 J[0][1] = X[2][0] - X[0][0];
646 J[0][2] = X[3][0] - X[0][0];
648 J[1][0] = X[1][1] - X[0][1];
649 J[1][1] = X[2][1] - X[0][1];
650 J[1][2] = X[3][1] - X[0][1];
652 J[2][0] = X[1][2] - X[0][2];
653 J[2][1] = X[2][2] - X[0][2];
654 J[2][2] = X[3][2] - X[0][2];
655 real64 const detJ = LvArray::tensorOps::determinant< 3 >( J );
656 return detJ * quadratureWeight( q );
659 template<
typename NUM_Q_POINTS >
664 real64 const (&X)[numNodes][3],
668 return calcJacobian( q, X, J );
673 template<
typename NUM_Q_POINTS >
679 real64 const (&X)[numNodes][3],
680 real64 (& gradN)[numNodes][3] )
683 gradN[0][0] = X[1][1]*( X[3][2] - X[2][2] ) - X[2][1]*( X[3][2] - X[1][2] ) + X[3][1]*( X[2][2] - X[1][2] );
684 gradN[0][1] = -X[1][0]*( X[3][2] - X[2][2] ) + X[2][0]*( X[3][2] - X[1][2] ) - X[3][0]*( X[2][2] - X[1][2] );
685 gradN[0][2] = X[1][0]*( X[3][1] - X[2][1] ) - X[2][0]*( X[3][1] - X[1][1] ) + X[3][0]*( X[2][1] - X[1][1] );
687 gradN[1][0] = -X[0][1]*( X[3][2] - X[2][2] ) + X[2][1]*( X[3][2] - X[0][2] ) - X[3][1]*( X[2][2] - X[0][2] );
688 gradN[1][1] = X[0][0]*( X[3][2] - X[2][2] ) - X[2][0]*( X[3][2] - X[0][2] ) + X[3][0]*( X[2][2] - X[0][2] );
689 gradN[1][2] = -X[0][0]*( X[3][1] - X[2][1] ) + X[2][0]*( X[3][1] - X[0][1] ) - X[3][0]*( X[2][1] - X[0][1] );
691 gradN[2][0] = X[0][1]*( X[3][2] - X[1][2] ) - X[1][1]*( X[3][2] - X[0][2] ) + X[3][1]*( X[1][2] - X[0][2] );
692 gradN[2][1] = -X[0][0]*( X[3][2] - X[1][2] ) + X[1][0]*( X[3][2] - X[0][2] ) - X[3][0]*( X[1][2] - X[0][2] );
693 gradN[2][2] = X[0][0]*( X[3][1] - X[1][1] ) - X[1][0]*( X[3][1] - X[0][1] ) + X[3][0]*( X[1][1] - X[0][1] );
695 gradN[3][0] = -X[0][1]*( X[2][2] - X[1][2] ) + X[1][1]*( X[2][2] - X[0][2] ) - X[2][1]*( X[1][2] - X[0][2] );
696 gradN[3][1] = X[0][0]*( X[2][2] - X[1][2] ) - X[1][0]*( X[2][2] - X[0][2] ) + X[2][0]*( X[1][2] - X[0][2] );
697 gradN[3][2] = -X[0][0]*( X[2][1] - X[1][1] ) + X[1][0]*( X[2][1] - X[0][1] ) - X[2][0]*( X[1][1] - X[0][1] );
699 real64 detJ = determinantJacobianTransformation( X );
700 real64 factor = 1.0 / ( detJ );
702 for(
int i = 0; i < numNodes; ++i )
704 for(
int j = 0; j < 3; ++j )
706 gradN[i][j] *= factor;
710 return detJ * weight * quadratureWeight( q );
713 template<
typename NUM_Q_POINTS >
718 real64 const (&X)[numNodes][3],
720 real64 ( & gradN )[numNodes][3] )
722 return calcGradN( q, X, gradN );
725 template<
typename NUM_Q_POINTS >
730 real64 const (&X)[numNodes][3],
731 real64 (& gradN)[numFaces][3] )
738 J[0][0] = X[1][0] - X[0][0];
739 J[0][1] = X[2][0] - X[0][0];
740 J[0][2] = X[3][0] - X[0][0];
742 J[1][0] = X[1][1] - X[0][1];
743 J[1][1] = X[2][1] - X[0][1];
744 J[1][2] = X[3][1] - X[0][1];
746 J[2][0] = X[1][2] - X[0][2];
747 J[2][1] = X[2][2] - X[0][2];
748 J[2][2] = X[3][2] - X[0][2];
750 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
752 real64 dNdXi[numFaces][3] = {{0}};
754 real64 const r = quadratureParentCoords0( q );
755 real64 const s = quadratureParentCoords1( q );
756 real64 const t = quadratureParentCoords2( q );
758 dNdXi[0][0] = ( 1 - 2 * r - s - t ) * t;
759 dNdXi[0][1] = -r * t;
760 dNdXi[0][2] = ( 1 - r - s - 2 * t ) * r;
762 dNdXi[1][0] = ( 1 - 2 * r - s - t ) * s;
763 dNdXi[1][1] = ( 1 - r - 2 * s - t ) * r;
764 dNdXi[1][2] = -r * s;
766 dNdXi[2][0] = -t * s;
767 dNdXi[2][1] = ( 1 - r - 2 * s - t ) * t;
768 dNdXi[2][2] = ( 1 - r - s - 2 * t ) * s;
774 for(
int fi=0; fi<numFaces; ++fi )
776 gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
777 gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
778 gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
781 return detJ * weight * quadratureWeight( q );
786 template<
typename NUM_Q_POINTS >
792 real64 const (&X)[numNodes][3] )
795 real64 detJ = determinantJacobianTransformation( X );
797 return detJ * weight * quadratureWeight( q );
804 template<
typename NUM_Q_POINTS >
#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 localIndex numSupportPoints
The number of support points per element.
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 calcGradFaceBubbleN(localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numFaces][3])
Calculate the bubble function derivatives wrt the physical coordinates.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcFaceBubbleN(real64 const (&pointCoord)[3], real64(&N)[numFaces])
Calculate face bubble functions values for each face at a given point in the parent space.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate face bubble functions values for each face at a quadrature point.
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.
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 localIndex getMaxSupportPoints()
Get the maximum number of support points.
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.
constexpr static int numSamplingPointsPerDirection
Number of sampling points per direction.
constexpr static localIndex numFaces
Number of faces in the 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 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 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature point.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints()
Get the number of support points.
static GEOS_HOST_DEVICE localIndex getNumSupportPoints(StackVariables 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.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN(real64 const (&pointCoord)[3], real64(&N)[numNodes])
Calculate shape functions values for each support point at a given point in the parent space.
constexpr static localIndex numSupportPoints
Number of support points in the element.
DO_NOT_DOCUMENT static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints()
Get the number of quadrature points.
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 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.
constexpr static int numSamplingPoints
Total number of sampling points.
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.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > * getImpl()
Get the device-compatible implementation type.
H1_Tetrahedron_Lagrange1_Gauss()
Constructor.
H1_Tetrahedron_Lagrange1_Gauss_impl< NUM_Q_POINTS > const * getImpl() const
Get the device-compatible implementation type.
virtual ~H1_Tetrahedron_Lagrange1_Gauss() override final=default
Destructor.
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.