20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_TRILINEARHEXAHEDRON_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_TRILINEARHEXAHEDRON_HPP_
84 template<
typename SubregionType >
159 real64 (& samplingPointCoord)[3] )
167 samplingPointCoord[0] = -1 + i0 * step;
168 samplingPointCoord[1] = -1 + i1 * step;
169 samplingPointCoord[2] = -1 + i2 * step;
223 return calcN( q, N );
380 return LvArray::tensorOps::invert< 3 >( J );
396 real64 const (&invJ)[3][3],
419 real64 const (&invJ)[3][3],
441 real64 const (&invJ)[3][3],
480 real64 const ( &invJ )[3][3],
489 constexpr
static real64 parentVolume = parentLength*parentLength*parentLength;
497 constexpr
static real64 quadratureFactor = 1.0 / 1.732050807568877293528;
518 template<
typename FUNC,
typename ... PARAMS >
520 static void supportLoop(
int const qa,
524 PARAMS &&... params );
543 template<
typename FUNC,
typename ... PARAMS >
545 H1_Hexahedron_Lagrange1_GaussLegendre2_impl::supportLoop(
int const qa,
549 PARAMS &&... params )
553 #define PARENT_GRADIENT_METHOD 2
554 #if PARENT_GRADIENT_METHOD == 1
558 real64 const quadratureCoords[3] = { -quadratureFactor + 1.154700538379252 * qa,
559 -quadratureFactor + 1.154700538379252 * qb,
560 -quadratureFactor + 1.154700538379252 * qc };
562 real64 const psi0[2] = { 0.5 - 0.5 * quadratureCoords[0],
563 0.5 + 0.5 * quadratureCoords[0] };
564 real64 const psi1[2] = { 0.5 - 0.5 * quadratureCoords[1],
565 0.5 + 0.5 * quadratureCoords[1] };
566 real64 const psi2[2] = { 0.5 - 0.5 * quadratureCoords[2],
567 0.5 + 0.5 * quadratureCoords[2] };
568 constexpr
real64 dpsi[2] = { -0.5, 0.5 };
569 #elif PARENT_GRADIENT_METHOD == 2
589 constexpr
static real64 psiProduct[3] = { 0.311004233964073108, 0.083333333333333333, 0.022329099369260226};
590 constexpr
static int dpsi[2] = { -1, 1 };
597 for(
int a=0; a<2; ++a )
599 #if PARENT_GRADIENT_METHOD == 2
600 int const qaa = ( a^qa );
602 for(
int b=0; b<2; ++b )
604 #if PARENT_GRADIENT_METHOD == 2
605 int const qbb = ( b^qb );
607 for(
int c=0; c<2; ++c )
609 #if PARENT_GRADIENT_METHOD == 2
610 int const qcc = ( c^qc );
613 #if PARENT_GRADIENT_METHOD == 1
614 real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2[c],
615 psi0[a] * dpsi[b] * psi2[c],
616 psi0[a] * psi1[b] * dpsi[c] };
617 #elif PARENT_GRADIENT_METHOD == 2
624 real64 const dNdXi[3] = { dpsi[a] * psiProduct[ qbb + qcc ],
625 dpsi[b] * psiProduct[ qaa + qcc ],
626 dpsi[c] * psiProduct[ qaa + qbb ] };
630 func( dNdXi, nodeIndex, std::forward< PARAMS >( params )... );
641 real64 const (&X)[numNodes][3],
644 for(
int i = 0; i < 3; ++i )
646 for(
int j = 0; j < 3; ++j )
654 real64 const detJ = LvArray::tensorOps::determinant< 3 >( J );
662 real64 const (&X)[numNodes][3],
673 real64 const (&X)[numNodes][3],
674 real64 (& gradN)[numNodes][3] )
683 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
687 return detJ * weight;
694 real64 const (&X)[numNodes][3],
696 real64 ( & gradN )[numNodes][3] )
705 real64 const (&X)[numNodes][3],
706 real64 (& gradN)[numFaces][3] )
716 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
728 gradN[fi][0] = dNdXi[fi][0] * J[0][0] + dNdXi[fi][1] * J[1][0] + dNdXi[fi][2] * J[2][0];
729 gradN[fi][1] = dNdXi[fi][0] * J[0][1] + dNdXi[fi][1] * J[1][1] + dNdXi[fi][2] * J[2][1];
730 gradN[fi][2] = dNdXi[fi][0] * J[0][2] + dNdXi[fi][1] * J[1][2] + dNdXi[fi][2] * J[2][2];
733 return detJ * weight;
738 #pragma GCC diagnostic push
739 #pragma GCC diagnostic ignored "-Wshadow"
749 real64 const (&X)[numNodes][3],
758 for(
int i = 0; i < 3; ++i )
760 for(
int j = 0; j < 3; ++j )
762 J[i][j] = J[i][j] + dNdXi[ j ] * Xnode[i];
788 real64 const ( &invJ )[3][3],
789 real64 (& gradN)[numNodes][3] )
793 real64 const (&invJ)[3][3],
805 gradN[nodeIndex][0] = dNdXi[0] * invJ[0][0] + dNdXi[1] * invJ[1][0] + dNdXi[2] * invJ[2][0];
806 gradN[nodeIndex][1] = dNdXi[0] * invJ[0][1] + dNdXi[1] * invJ[1][1] + dNdXi[2] * invJ[2][1];
807 gradN[nodeIndex][2] = dNdXi[0] * invJ[0][2] + dNdXi[1] * invJ[1][2] + dNdXi[2] * invJ[2][2];
819 real64 const (&X)[numNodes][3] )
828 return LvArray::tensorOps::determinant< 3 >( J );
836 real64 const (&invJ)[3][3],
837 real64 const (&var)[numNodes][3],
845 real64 const (&invJ)[3][3],
850 real64 gradN[3] = {0, 0, 0};
851 for(
int i = 0; i < 3; ++i )
853 for(
int j = 0; j < 3; ++j )
855 gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
859 grad[0] = grad[0] + gradN[0] * var[ nodeIndex ][0];
860 grad[1] = grad[1] + gradN[1] * var[ nodeIndex ][1];
861 grad[2] = grad[2] + gradN[2] * var[ nodeIndex ][2];
862 grad[3] = grad[3] + gradN[2] * var[ nodeIndex ][1] + gradN[1] * var[ nodeIndex ][2];
863 grad[4] = grad[4] + gradN[2] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][2];
864 grad[5] = grad[5] + gradN[1] * var[ nodeIndex ][0] + gradN[0] * var[ nodeIndex ][1];
865 }, invJ, var, grad );
871 real64 const (&invJ)[3][3],
873 real64 (& R)[numNodes][3] )
878 supportLoop( qa, qb, qc,
880 (
real64 const (&dNdXi)[3],
882 real64 const (&invJ)[3][3],
887 real64 gradN[3] = {0, 0, 0};
888 for(
int i = 0; i < 3; ++i )
890 for(
int j = 0; j < 3; ++j )
892 gradN[i] = gradN[i] + dNdXi[ j ] * invJ[j][i];
895 R[ nodeIndex ][ 0 ] = R[ nodeIndex ][ 0 ] - var[ 0 ] * gradN[ 0 ] - var[ 5 ] * gradN[ 1 ] - var[ 4 ] * gradN[ 2 ];
896 R[ nodeIndex ][ 1 ] = R[ nodeIndex ][ 1 ] - var[ 5 ] * gradN[ 0 ] - var[ 1 ] * gradN[ 1 ] - var[ 3 ] * gradN[ 2 ];
897 R[ nodeIndex ][ 2 ] = R[ nodeIndex ][ 2 ] - var[ 4 ] * gradN[ 0 ] - var[ 3 ] * gradN[ 1 ] - var[ 2 ] * gradN[ 2 ];
906 real64 const (&invJ)[3][3],
907 real64 const (&var)[numNodes][3],
915 real64 const (&invJ)[3][3],
919 for(
int i = 0; i < 3; ++i )
922 for(
int j = 0; j < 3; ++j )
924 gradN = gradN + dNdXi[ j ] * invJ[j][i];
926 for(
int k = 0; k < 3; ++k )
928 grad[k][i] = grad[k][i] + gradN * var[ nodeIndex ][k];
931 }, invJ, var, grad );
937 #pragma GCC diagnostic pop
980 #undef PARENT_GRADIENT_METHOD
#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_RESTRICT
preprocessor variable for the C99 restrict keyword for use with pointers
#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 void symmetricGradient(int const q, real64 const (&invJ)[3][3], real64 const (&var)[numNodes][3], real64(&grad)[6])
Calculate the symmetric gradient of a vector valued support field at a quadrature point using the sto...
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()
Get the number of quadrature 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.
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 jacobianTransformation(int const qa, int const qb, int const qc, 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], StackVariables const &stack)
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], StackVariables const &stack, real64(&J)[3][3])
Calculate the Jacobian matrix at a quadrature 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.
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.
static GEOS_HOST_DEVICE void plusGradNajAij(int const q, real64 const (&invJ)[3][3], real64 const (&var)[6], real64(&R)[numNodes][3])
Inner product of all basis function gradients and a rank-2 symmetric tensor evaluated at a quadrature...
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 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 gradient(int const q, real64 const (&invJ)[3][3], real64 const (&var)[numNodes][3], real64(&grad)[3][3])
Calculate the gradient of a vector valued support field at a point using the stored basis function gr...
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
static GEOS_HOST_DEVICE void applyTransformationToParentGradients(int const qa, int const qb, int const qc, real64 const (&invJ)[3][3], real64(&gradN)[numNodes][3])
Apply a Jacobian transformation matrix from the parent space to the physical space on the parent shap...
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.
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.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
static GEOS_HOST_DEVICE localIndex getMaxSupportPoints()
Get the maximum number of support points.
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 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 bubble function derivatives wrt the physical coordinates.
H1_Hexahedron_Lagrange1_GaussLegendre2_impl const * getImpl() const
Get the device-compatible implementation type.
H1_Hexahedron_Lagrange1_GaussLegendre2_impl * getImpl()
Get the device-compatible implementation type.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 parentSupportCoord(const localIndex supportPointIndex)
Calculate the parent coordinates for the xi0 direction, given the linear index of a support point.
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.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void value(real64 const (&coords)[3], real64(&N)[numSupportPoints])
The value of the basis function for a support point evaluated at a point along the axes.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE void multiIndex(const int linearIndex, int &i0, int &i1, int &i2)
Calculate the Cartesian/TensorProduct index given the linear index of a support point.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void gradientFaceBubble(real64 const (&coords)[3], real64(&dNdXi)[numSupportFaces][3])
The value of the bubble basis function derivatives for a support face evaluated at a point along the ...
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE int linearIndex(const int i, const int j, const int k)
Calculates the linear index for support/quadrature points from ijk coordinates.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void valueFaceBubble(real64 const (&coords)[3], real64(&N)[numSupportFaces])
The value of the bubble basis function for a support face evaluated at a point along the axes.