20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_PKPYRAMIDBCD_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_PKPYRAMIDBCD_HPP_
24 #include "denseLinearAlgebra/interfaces/blaslapack/BlasLapackLA.hpp"
37 constexpr
int Pk_Pyramid_BCD_NumNodes = ( ORDER + 1 ) * ( ORDER + 2 ) * ( 2 * ORDER + 3 ) / 6;
50 Pk_Pyramid_BCD_NumNodes< ORDER > >
56 Pk_Pyramid_BCD_NumNodes< ORDER > >;
62 template<
typename SubregionType >
75 static constexpr
int order = ORDER;
158 for(
int l = 0; l < qc; ++l )
160 int n =
order + 1 - l;
163 int n_k =
order + 1 - qc;
164 index += qa + qb * n_k;
179 if constexpr (ORDER == 1)
181 coords[0][0] = -1.0; coords[0][1] = 1.0; coords[0][2] = 0.0;
182 coords[1][0] = -1.0; coords[1][1] = -1.0; coords[1][2] = 0.0;
183 coords[2][0] = 1.0; coords[2][1] = -1.0; coords[2][2] = 0.0;
184 coords[3][0] = 1.0; coords[3][1] = 1.0; coords[3][2] = 0.0;
185 coords[4][0] = 0.0; coords[4][1] = 0.0; coords[4][2] = 1.0;
187 else if constexpr (ORDER == 2)
199 template<
typename FUNC >
209 localIndex maxIj = LvArray::math::max( i, j );
210 for(
localIndex k = 0; k <= ORDER - maxIj; ++k )
222 template<
typename FUNC >
231 localIndex maxIj = LvArray::math::max( i, j );
232 for(
localIndex k = 0; k <= ORDER - maxIj; ++k )
251 for(
int i = 0; i < 4; i++ )
253 for(
int j = 0; j < 3; j++ )
255 m[ i ][ j ] = X[ i + 1 ][ j ] - X[ 0 ][ j ];
258 return LvArray::math::abs( LvArray::tensorOps::determinant< 3 >( m ) );
273 int i1 = ( face + 1 ) % 4;
274 int i2 = ( face + 2 ) % 4;
275 int i3 = ( face + 3 ) % 4;
277 real64 ab[3] = { X[ i2 ][ 0 ] - X[ i1 ][ 0 ],
278 X[ i2 ][ 1 ] - X[ i1 ][ 1 ],
279 X[ i2 ][ 2 ] - X[ i1 ][ 2 ] };
280 real64 ac[3] = { X[ i3 ][ 0 ] - X[ i1 ][ 0 ],
281 X[ i3 ][ 1 ] - X[ i1 ][ 1 ],
282 X[ i3 ][ 2 ] - X[ i1 ][ 2 ] };
283 real64 term1 = ab[1] * ac[2] - ab[2] * ac[1];
284 real64 term2 = ab[2] * ac[0] - ab[0] * ac[2];
285 real64 term3 = ab[0] * ac[1] - ab[1] * ac[0];
286 return LvArray::math::sqrt( ( term1 * term1 + term2 * term2 + term3 * term3 ) );
305 return calcN( q, N );
330 real64 val = 0.5 * (alpha - beta + (alpha + beta + 2) * x);
336 real64 P_prev1 = 0.5 * (alpha - beta + (alpha + beta + 2) * x);
340 for(
int k = 2; k < n+1; ++k )
342 real64 a_k = 2 * k * (k + alpha + beta) * (2 * k + alpha + beta - 2);
343 real64 b_k = (2 * k + alpha + beta - 1) * (alpha * alpha - beta * beta);
344 real64 c_k = (2 * k + alpha + beta - 1) * (2 * k + alpha + beta) * (2 * k + alpha + beta - 2);
345 real64 d_k = 2 * (k + alpha - 1) * (k + beta - 1) * (2 * k + alpha + beta);
347 P_current = ((b_k + c_k * x) * P_prev1 - d_k * P_prev2) / a_k;
376 real64 coeff = 0.5 * (n + alpha + beta + 1.0);
379 return coeff * jacobi_nm1;
398 real64 const epsilon = 1e-10;
399 if( LvArray::math::abs( z-1.0 ) < epsilon )
401 if( i == 0 && j == 0 )
403 return (k+2)*(k+1)/2.0;
411 real64 xi = x / (1.0 - z);
412 real64 eta = y / (1.0 - z);
413 real64 chi = 2.0 * z - 1.0;
414 localIndex max_ij = LvArray::math::max( i, j );
419 return P_i * P_j * std::pow( 1.0 - z, max_ij ) * P_k;
439 real64 const epsilon = 1e-10;
440 if( LvArray::math::abs( z-1.0 ) < epsilon )
449 real64 xi = x / (1.0 - z);
450 real64 eta = y / (1.0 - z);
451 real64 chi = 2.0 * z - 1.0;
460 real64 f = pow( 1.0 - z, m );
462 gradPsiX[0] = (1.0 / (1.0 - z)) * dP_i * P_j * f * P_k;
465 gradPsiX[1] = (1.0 / (1.0 - z)) * dP_j * P_i * f * P_k;
468 real64 dxi_dz = x / std::pow( 1.0 - z, 2 );
469 real64 deta_dz = y / std::pow( 1.0 - z, 2 );
471 real64 df_dz = -m * std::pow( 1.0 - z, m - 1 );
473 dP_i * dxi_dz * P_j * f * P_k +
474 P_i * dP_j * deta_dz * f * P_k +
475 P_i * P_j * df_dz * P_k +
476 P_i * P_j * f * dP_k * dchi_dz;
500 PsiX[count] =
calcModal( p, r, s, coords[j] );
501 VDM[count][j] = PsiX[count];
508 BlasLapackLA::matrixInverse( VDM, VDM_inv );
534 PsiX[count] =
calcModal( p, r, s, coords[q] );
544 N[i] += VDM_inv[i][j] * PsiX[j];
576 N[i] += VDM_inv[i][j] * PsiX[j];
616 gradN[i][0] += VDM_inv[i][j] * gradModal[j][0];
617 gradN[i][1] += VDM_inv[i][j] * gradModal[j][1];
618 gradN[i][2] += VDM_inv[i][j] * gradModal[j][2];
653 gradN[i][0] += VDM_inv[i][j] * gradModal[j][0];
654 gradN[i][1] += VDM_inv[i][j] * gradModal[j][1];
655 gradN[i][2] += VDM_inv[i][j] * gradModal[j][2];
668 template<
typename FUNC >
675 constexpr
real64 GLeQuadraturePoints[3] = { -0.7745966692, 0.0, 0.7745966692 };
676 constexpr
real64 GLeQuadratureWeights[3] = { 0.5555555556, 0.8888888889, 0.5555555556 };
678 constexpr
real64 GJQuadraturePoints[3] = { 0.07299, 0.34700, 0.70500 };
679 constexpr
real64 GJQuadratureWeights[3] = { 0.15714, 0.14625, 0.02995 };
689 for(
int i = 0; i < 3; ++i )
691 for(
int j = 0; j < 3; ++j )
693 for(
int k = 0; k < 3; ++k )
696 real64 xi = GLeQuadraturePoints[i];
697 real64 eta = GLeQuadraturePoints[j];
698 real64 chi = GJQuadraturePoints[k];
700 real64 weight = GLeQuadratureWeights[i] * GLeQuadratureWeights[j] * GJQuadratureWeights[k];
704 real64 X[3] = {x_i, y_j, z_k};
706 calcN( X, VDM_inv, N );
708 val += weight * N[a] * N[b];
725 template<
typename FUNC >
731 constexpr
real64 GLeQuadraturePoints[3] = { -0.7745966692, 0.0, 0.7745966692 };
732 constexpr
real64 GLeQuadratureWeights[3] = { 0.5555555556, 0.8888888889, 0.5555555556 };
734 constexpr
real64 GJQuadraturePoints[3] = { 0.07299, 0.34700, 0.70500 };
735 constexpr
real64 GJQuadratureWeights[3] = { 0.15714, 0.14625, 0.02995 };
746 for(
int i = 0; i < 3; ++i )
748 for(
int j = 0; j < 3; ++j )
750 for(
int k = 0; k < 3; ++k )
752 real64 xi = GLeQuadraturePoints[i];
753 real64 eta = GLeQuadraturePoints[j];
754 real64 chi = GJQuadraturePoints[k];
756 real64 weight = GLeQuadratureWeights[i] * GLeQuadratureWeights[j] * GJQuadratureWeights[k];
760 real64 X[3] = {x_i, y_j, z_k};
763 real64 dot = gradN[a][0] * gradN[b][0]
764 + gradN[a][1] * gradN[b][1]
765 + gradN[a][2] * gradN[b][2];
778 template<
int ORDER >
803 #pragma nv_diag_suppress 20012
817 return static_cast< ImplType *
>(
this);
826 return static_cast< ImplType const *
>(
this);
#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.
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 maxSupportPoints
The maximum number of support points per element.
constexpr static localIndex numNodes
The number of nodes per element.
Base class for FEM element implementations.
static GEOS_FORCE_INLINE array2d< real64 > computeVanderMondeMatrixInverse()
Compute the inverse of the VanDerMonde matrix.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumSupportPoints(StackVariables const &stack)
Get the number of support points.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 faceJacobianDeterminant(localIndex face, real64 const (&X)[5][3])
Calculate the determinant of the jacobian on the face opposite to the given vertex.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void calcGradModal(int i, int j, int k, real64 const (&X)[3], real64(&gradPsiX)[3])
Calculate modal base functions derivative values for a modal point (i,j,k) at a point X.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumModes(StackVariables const &stack)
Get the number of modal points.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void calcGradN(localIndex q, arrayView2d< real64 const > VDM_inv, real64(&gradN)[numNodes][3])
Evaluate shape functions of a linear pyramid (5-node) at a quadrature point.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 EvaluateJacobiPolynomial(localIndex const n, real64 const alpha, real64 const beta, real64 const x)
Helper to compute Jacobi Polynomials.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumSupportPoints()
Get the number of support points.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getMaxSupportPoints()
Get the maximum number of support points.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void calcN(localIndex q, arrayView2d< real64 const > VDM_inv, real64(&N)[numNodes])
Evaluate shape functions of a linear pyramid (5-node) at a quadrature point.
static constexpr GEOS_FORCE_INLINE void generateIndexesHost(FUNC &&func)
Generate the indexes for the modal shape functions.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void calcGradN(real64 const (&X)[3], arrayView2d< real64 const > VDM_inv, real64(&gradN)[numNodes][3])
Evaluate shape functions of a linear pyramid (5-node) at a given point.
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE localIndex linearIndex3DVal(const localIndex qa, localIndex const qb, localIndex const qc)
The linear index associated to the given one-dimensional indices in the three directions.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void generatePointsCoordinates(real64(&coords)[numNodes][3])
Generate the coordinates of the support points.
static constexpr int order
The order of the finite element.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 EvaluateJacobiPolynomialDerivative(localIndex const n, real64 const alpha, real64 const beta, real64 const x)
Evaluate the derivative of the Jacobi polynomial at a point x.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN(localIndex const q, StackVariables const &stack, real64(&N)[numNodes])
Calculate shape functions values for each support point at a quadrature/nodal point.
static GEOS_HOST_DEVICE localIndex getNumQuadraturePoints()
Get the number of quadrature points.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumQuadraturePoints(StackVariables const &stack)
Get the number of quadrature points.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 calcModal(localIndex const i, localIndex const j, localIndex const k, real64 const (&X)[3])
Calculate modal base functions values for a modal point (i,j,k) at a point X.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void computeMassTerm(arrayView2d< real64 const > VDM_inv, FUNC &&func)
Computes the mass term Mij in the mass matrix.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void generateIndexes(FUNC &&func)
Generate the indexes for the modal shape functions.
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void calcN(real64 const (&X)[3], arrayView2d< real64 const > VDM_inv, real64(&N)[numNodes])
Evaluate shape functions of a linear pyramid (5-node) at a quadrature point.
constexpr static localIndex numModes
The number of modal points per element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 jacobianDeterminant(real64 const (&X)[5][3])
Returns the determinant of the Jacobian of the element.
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void computeStiffnessTerm(arrayView2d< real64 const > VDM_inv, FUNC &&func)
Computes the stifness term Kij in the mass matrix.
typename Base::template MeshData< SubregionType > MeshData
Mesh data structure for the element.
ImplType const * getImpl() const
Get the device-compatible implementation type.
ImplType * getImpl()
Get the device-compatible implementation type.
Array< T, 2, PERMUTATION > array2d
Alias for 2D array.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Kernel variables (dof numbers, jacobian and residual) located on the stack.