20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1PYRAMIDLAGRANGE1GAUSS5_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_H1PYRAMIDLAGRANGE1GAUSS5_HPP_
135 real64 (& samplingPointCoord)[3] )
189 GEOS_ERROR(
"Unsupported bubble functions for pyramid elements" );
205 GEOS_ERROR(
"Unsupported bubble functions for pyramid elements" );
291 jacobianTransformation( q, X, J );
292 return LvArray::tensorOps::invert< 3 >( J );
297 constexpr
static real64 weight = 81.0 / 100.0;
300 constexpr
static real64 weightDelta = 125.0 / 27.0 - weight;
305 constexpr
static real64 quadratureCrossSectionCoord = 0.584237394672177;
309 constexpr
static real64 quadratureLongitudinalCoordNeg = -2.0 / 3.0;
313 constexpr
static real64 quadratureLongitudinalCoordDelta = 16.0 / 15.0;
322 template<
typename T >
325 constexpr
static T linearMap( T
const i, T
const j )
340 return -1.0 + 2.0 * ( a & 1 ) + 0.25 * ( a & 4 );
353 return -1.0 + ( a & 2 ) + 0.25 * ( a & 4 );
366 return -1.0 + 0.5 * ( a & 4 );
379 return parentCoords0( q ) * quadratureCrossSectionCoord;
392 return parentCoords1( q ) * quadratureCrossSectionCoord;
405 return quadratureLongitudinalCoordNeg + 0.5 * ( 1 + parentCoords2( q ) ) * quadratureLongitudinalCoordDelta;
417 return weight + 0.5 * ( 1 + parentCoords2( q ) ) * weightDelta;
427 static void jacobianTransformation(
int const q,
442 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
443 real64 const ( &invJ )[3][3],
453 H1_Pyramid_Lagrange1_Gauss5::
454 jacobianTransformation(
int const q,
455 real64 const (&X)[numNodes][3],
458 real64 const quadratureCoords[3] = { quadratureParentCoords0( q ),
459 quadratureParentCoords1( q ),
460 quadratureParentCoords2( q ) };
462 real64 const psi0[2] = { 0.5 - 0.5*quadratureCoords[0],
463 0.5 + 0.5*quadratureCoords[0] };
464 real64 const psi1[2] = { 0.5 - 0.5*quadratureCoords[1],
465 0.5 + 0.5*quadratureCoords[1] };
466 real64 const psi2 = 0.5 - 0.5*quadratureCoords[2];
467 constexpr
real64 dpsi[2] = { -0.5, 0.5 };
474 real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2,
475 psi0[a] * dpsi[b] * psi2,
476 psi0[a] * psi1[b] * dpsi[0] };
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];
489 for(
int i = 0; i < 3; ++i )
491 J[i][2] = J[i][2] + dpsi[1] * X[4][i];
500 H1_Pyramid_Lagrange1_Gauss5::
501 applyJacobianTransformationToShapeFunctionsDerivatives(
int const q,
502 real64 const ( &invJ )[3][3],
503 real64 (& gradN)[numNodes][3] )
505 real64 const quadratureCoords[3] = { quadratureParentCoords0( q ),
506 quadratureParentCoords1( q ),
507 quadratureParentCoords2( q ) };
509 real64 const psi0[2] = { 0.5*( 1.0 - quadratureCoords[0] ),
510 0.5*( 1.0 + quadratureCoords[0] ) };
511 real64 const psi1[2] = { 0.5*( 1.0 - quadratureCoords[1] ),
512 0.5*( 1.0 + quadratureCoords[1] ) };
513 real64 const psi2 = 0.5*( 1.0 - quadratureCoords[2]);
514 constexpr
real64 dpsi[2] = { -0.5, 0.5 };
521 real64 const dNdXi[3] = { dpsi[a] * psi1[b] * psi2,
522 psi0[a] * dpsi[b] * psi2,
523 psi0[a] * psi1[b] * dpsi[0] };
524 localIndex const nodeIndex = linearMap( a, b );
525 for(
int i = 0; i < 3; ++i )
527 gradN[nodeIndex][i] = 0.0;
528 for(
int j = 0; j < 3; ++j )
530 gradN[nodeIndex][i] = gradN[nodeIndex][i] + dNdXi[ j ] * invJ[j][i];
537 for(
int i = 0; i < 3; ++i )
539 gradN[4][i] = dpsi[1] * invJ[2][i];
550 real64 ( & N )[numNodes] )
552 N[0] = 0.125*( 1.0 - pointCoord[0] ) * ( 1.0 - pointCoord[1] ) * ( 1.0 - pointCoord[2] );
553 N[1] = 0.125*( 1.0 + pointCoord[0] ) * ( 1.0 - pointCoord[1] ) * ( 1.0 - pointCoord[2] );
554 N[2] = 0.125*( 1.0 - pointCoord[0] ) * ( 1.0 + pointCoord[1] ) * ( 1.0 - pointCoord[2] );
555 N[3] = 0.125*( 1.0 + pointCoord[0] ) * ( 1.0 + pointCoord[1] ) * ( 1.0 - pointCoord[2] );
556 N[4] = 0.5*( 1.0 + pointCoord[2] );
564 real64 ( & N )[numNodes] )
566 real64 const xi[3] = { quadratureParentCoords0( q ),
567 quadratureParentCoords1( q ),
568 quadratureParentCoords2( q ) };
578 real64 ( & N )[numNodes] )
580 return calcN( q, N );
588 real64 const (&X)[numNodes][3],
589 real64 (& gradN)[numNodes][3] )
593 jacobianTransformation( q, X, J );
595 real64 const detJ = LvArray::tensorOps::invert< 3 >( J );
597 applyJacobianTransformationToShapeFunctionsDerivatives( q, J, gradN );
599 return detJ * quadratureWeight( q );
606 real64 const (&X)[numNodes][3],
608 real64 ( & gradN )[numNodes][3] )
617 real64 const (&X)[numNodes][3],
618 real64 (& gradN)[numFaces][3] )
621 GEOS_ERROR(
"Unsupported bubble functions for pyramid elements" );
632 real64 const (&X)[numNodes][3] )
636 jacobianTransformation( q, X, J );
638 return LvArray::tensorOps::determinant< 3 >( J ) * 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.
#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.
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const override
Virtual getter for the number of support points per element.
constexpr static localIndex numFaces
The number of faces/support 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.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
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, 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 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.
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 calcN(real64 const (&pointCoord)[3], real64(&N)[numNodes])
Calculate shape functions values at a single 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.
constexpr static localIndex numNodes
The number of nodes/support points per element.
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.
constexpr static localIndex numQuadraturePoints
The number of quadrature 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...
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void getSamplingPointCoordInParentSpace(int const &linearIndex, real64(&samplingPointCoord)[3])
Get the Sampling Point Coord In the Parent Space.
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.
static GEOS_HOST_DEVICE void calcFaceBubbleN(localIndex const q, real64(&N)[numFaces])
Calculate shape functions values for each support face at a quadrature point.
static GEOS_HOST_DEVICE real64 transformedQuadratureWeight(localIndex const q, real64 const (&X)[numNodes][3])
Calculate the integration weights for a quadrature point.
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.