GEOS
Pk_Pyramid_BCD.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2018-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_PKPYRAMIDBCD_HPP_
21 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_PKPYRAMIDBCD_HPP_
22 
23 #include "FiniteElementBase.hpp"
24 #include "denseLinearAlgebra/interfaces/blaslapack/BlasLapackLA.hpp"
25 #include <utility>
26 
27 
28 
29 namespace geos
30 {
31 namespace finiteElement
32 {
33 
34 namespace
35 {
36 template< int ORDER >
37 constexpr int Pk_Pyramid_BCD_NumNodes = ( ORDER + 1 ) * ( ORDER + 2 ) * ( 2 * ORDER + 3 ) / 6;
38 }
39 
40 
47 template< int ORDER >
48 class Pk_Pyramid_BCD_impl : public FiniteElementBase_impl< Pk_Pyramid_BCD_NumNodes< ORDER >,
49  5,
50  Pk_Pyramid_BCD_NumNodes< ORDER > >
51 {
52 public:
55  5,
56  Pk_Pyramid_BCD_NumNodes< ORDER > >;
57 
59  using StackVariables = typename Base::StackVariables;
60 
62  template< typename SubregionType >
63  using MeshData = typename Base::template MeshData< SubregionType >;
64 
66  using Base::numNodes;
67 
70 
73 
75  static constexpr int order = ORDER;
76 
78  constexpr static localIndex numModes = numNodes;
79 
83  {
84  return numQuadraturePoints;
85  }
86 
95  {
96  GEOS_UNUSED_VAR( stack );
97  return numQuadraturePoints;
98  }
99 
104  {
105  return numNodes;
106  }
107 
112  {
113  return maxSupportPoints;
114  }
115 
124  {
125  GEOS_UNUSED_VAR( stack );
126  return numNodes;
127  }
128 
129 
130 
138  static localIndex getNumModes( StackVariables const & stack )
139  {
140  GEOS_UNUSED_VAR( stack );
141  return numNodes;
142  }
143 
144 
154  constexpr static localIndex linearIndex3DVal( const localIndex qa, localIndex const qb, localIndex const qc )
155  {
156  localIndex index = 0;
157 
158  for( int l = 0; l < qc; ++l )
159  {
160  int n = order + 1 - l;
161  index += n * n;
162  }
163  int n_k = order + 1 - qc;
164  index += qa + qb * n_k;
165 
166  return index;
167  }
168 
169 
176  static constexpr void generatePointsCoordinates( real64 (& coords)[numNodes][3] )
177  {
178  // Generate the coordinates of the support points based on the order
179  if constexpr (ORDER == 1)
180  {
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;
186  }
187  else if constexpr (ORDER == 2)
188  {
189  //TODO
190  }
191 
192  }
193 
194 
199  template< typename FUNC >
202  static constexpr void generateIndexes( FUNC && func )
203  {
204 
205  for( localIndex i = 0; i <= ORDER; ++i )
206  {
207  for( localIndex j = 0; j <= ORDER; ++j )
208  {
209  localIndex maxIj = LvArray::math::max( i, j );
210  for( localIndex k = 0; k <= ORDER - maxIj; ++k )
211  {
212  func( i, j, k );
213  }
214  }
215  }
216  }
217 
222  template< typename FUNC >
224  static constexpr void generateIndexesHost( FUNC && func )
225  {
226 
227  for( localIndex i = 0; i <= ORDER; ++i )
228  {
229  for( localIndex j = 0; j <= ORDER; ++j )
230  {
231  localIndex maxIj = LvArray::math::max( i, j );
232  for( localIndex k = 0; k <= ORDER - maxIj; ++k )
233  {
234  func( i, j, k );
235  }
236  }
237  }
238  }
239 
247  static real64 jacobianDeterminant( real64 const (&X)[5][3] )
248 
249  {
250  real64 m[3][3] = {};
251  for( int i = 0; i < 4; i++ )
252  {
253  for( int j = 0; j < 3; j++ )
254  {
255  m[ i ][ j ] = X[ i + 1 ][ j ] - X[ 0 ][ j ];
256  }
257  }
258  return LvArray::math::abs( LvArray::tensorOps::determinant< 3 >( m ) );
259  }
260 
261 
271  real64 const (&X)[5][3] )
272  {
273  int i1 = ( face + 1 ) % 4;
274  int i2 = ( face + 2 ) % 4;
275  int i3 = ( face + 3 ) % 4;
276 
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 ) );
287 
288  }
289 
300  static void calcN( localIndex const q,
301  StackVariables const & stack,
302  real64 ( & N )[numNodes] )
303  {
304  GEOS_UNUSED_VAR( stack );
305  return calcN( q, N );
306  }
307 
308 
319  static constexpr real64 EvaluateJacobiPolynomial( localIndex const n, real64 const alpha, real64 const beta, real64 const x )
320  {
321 
322  if( n == 0 )
323  {
324  real64 val = 1.0;
325  return val;
326  }
327 
328  if( n == 1 )
329  {
330  real64 val = 0.5 * (alpha - beta + (alpha + beta + 2) * x);
331  return val;
332  }
333 
334  // Initial conditions for the recurrence relation
335  real64 P_prev2 = 1.0; // P_0
336  real64 P_prev1 = 0.5 * (alpha - beta + (alpha + beta + 2) * x); // P_1
337  real64 P_current = 0.0;
338 
339  // Recurrence relation to compute the Jacobi polynomial
340  for( int k = 2; k < n+1; ++k )
341  {
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);
346 
347  P_current = ((b_k + c_k * x) * P_prev1 - d_k * P_prev2) / a_k;
348 
349  // Mise à jour pour l'itération suivante
350  P_prev2 = P_prev1;
351  P_prev1 = P_current;
352  }
353 
354  return P_current;
355  }
356 
367  static constexpr real64 EvaluateJacobiPolynomialDerivative( localIndex const n, real64 const alpha, real64 const beta, real64 const x )
368  {
369  // Particular case for n = 0
370  if( n == 0 )
371  {
372  return 0.0;
373  }
374 
375  // dP_n^{(alpha,beta)}(x)/dx = (n+alpha+beta+1)/2 * P_{n-1}^{(alpha+1,beta+1)}(x)
376  real64 coeff = 0.5 * (n + alpha + beta + 1.0);
377  real64 jacobi_nm1 = EvaluateJacobiPolynomial( n-1, alpha+1.0, beta+1.0, x );
378 
379  return coeff * jacobi_nm1;
380  }
381 
391  static constexpr real64 calcModal( localIndex const i, localIndex const j, localIndex const k, real64 const (&X)[3] )
392  {
393 
394  real64 const x = X[0];
395  real64 const y = X[1];
396  real64 const z = X[2];
397 
398  real64 const epsilon = 1e-10; // Small value to avoid compairing floating point numbers directly
399  if( LvArray::math::abs( z-1.0 ) < epsilon )
400  {
401  if( i == 0 && j == 0 )
402  {
403  return (k+2)*(k+1)/2.0;
404  }
405  else
406  {
407  return 0.0;
408  }
409  }
410 
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 );
415  real64 P_i = EvaluateJacobiPolynomial( i, 0.0, 0.0, xi );
416  real64 P_j = EvaluateJacobiPolynomial( j, 0.0, 0.0, eta );
417  real64 P_k = EvaluateJacobiPolynomial( k, 2.0 * max_ij+2.0, 0.0, chi );
418 
419  return P_i * P_j * std::pow( 1.0 - z, max_ij ) * P_k;
420 
421  }
422 
432  static constexpr void calcGradModal( int i, int j, int k, real64 const (&X)[3],
433  real64 (& gradPsiX)[3] )
434  {
435  real64 x = X[0];
436  real64 y = X[1];
437  real64 z = X[2];
438 
439  real64 const epsilon = 1e-10; // Small value to avoid compairing floating point numbers directly
440  if( LvArray::math::abs( z-1.0 ) < epsilon )
441  {
442  gradPsiX[0] = 0.0;
443  gradPsiX[1] = 0.0;
444  gradPsiX[2] = 0.0;
445  return;
446  }
447 
448 
449  real64 xi = x / (1.0 - z);
450  real64 eta = y / (1.0 - z);
451  real64 chi = 2.0 * z - 1.0;
452  localIndex m = LvArray::math::max( i, j );
453  real64 P_i = EvaluateJacobiPolynomial( i, 0.0, 0.0, xi );
454  real64 P_j = EvaluateJacobiPolynomial( j, 0.0, 0.0, eta );
455  real64 P_k = EvaluateJacobiPolynomial( k, 2.0 * m + 2.0, 0.0, chi );
456  real64 dP_i = EvaluateJacobiPolynomialDerivative( i, 0.0, 0.0, xi );
457  real64 dP_j = EvaluateJacobiPolynomialDerivative( j, 0.0, 0.0, eta );
458  real64 dP_k = EvaluateJacobiPolynomialDerivative( k, 2.0 * m + 2.0, 0.0, chi );
459 
460  real64 f = pow( 1.0 - z, m );
461  // Gradient en x
462  gradPsiX[0] = (1.0 / (1.0 - z)) * dP_i * P_j * f * P_k;
463 
464  // Gradient en y
465  gradPsiX[1] = (1.0 / (1.0 - z)) * dP_j * P_i * f * P_k;
466 
467  // Gradient en z
468  real64 dxi_dz = x / std::pow( 1.0 - z, 2 );
469  real64 deta_dz = y / std::pow( 1.0 - z, 2 );
470  real64 dchi_dz = 2.0;
471  real64 df_dz = -m * std::pow( 1.0 - z, m - 1 );
472  gradPsiX[2] =
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;
477 
478 
479  }
480 
487  {
488  array2d< real64 > VDM;
489  VDM.resize( numNodes, numNodes );
490  real64 PsiX[numNodes] = {};
491  real64 coords[numNodes][3] = {};
492 
493  generatePointsCoordinates( coords );
494 
495  for( int j = 0; j < numNodes; ++j )
496  {
497  localIndex count = 0;
498  generateIndexesHost( [&]( localIndex const p, localIndex const r, localIndex const s )
499  {
500  PsiX[count] = calcModal( p, r, s, coords[j] );
501  VDM[count][j] = PsiX[count];
502  ++count;
503  } );
504  }
505  array2d< real64 > VDM_inv;
506  VDM_inv.resize( numNodes, numNodes );
507  // Inversion of VanDerMonde matrix
508  BlasLapackLA::matrixInverse( VDM, VDM_inv );
509  return VDM_inv;
510  }
511 
520  static constexpr void calcN( localIndex q, arrayView2d< real64 const > VDM_inv, real64 (& N)[numNodes] )
521  {
522 
523  real64 PsiX[numNodes] = {0.0};
524  real64 coords[numNodes][3] = {{0.0}};
525 
526  //Generate the coordinates of the support points
527  generatePointsCoordinates( coords );
528 
529 
530  localIndex count = 0;
531  //Compute the modal functions at the quadrature point
532  generateIndexes( [&]( localIndex const p, localIndex const r, localIndex const s )
533  {
534  PsiX[count] = calcModal( p, r, s, coords[q] );
535  ++count;
536  } );
537 
538  //Compute the nodal functions at the quadrature point
539  for( int i = 0; i < numNodes; ++i )
540  {
541  N[i] = 0.0;
542  for( int j = 0; j < numNodes; ++j )
543  {
544  N[i] += VDM_inv[i][j] * PsiX[j];
545  }
546  }
547 
548  }
557  static constexpr void calcN( real64 const (&X)[3], arrayView2d< real64 const > VDM_inv, real64 (& N)[numNodes] )
558  {
559 
560  real64 PsiX[numNodes] = {};
561 
562  localIndex count = 0;
563  //Compute the modal functions at the given point
564  generateIndexes( [&]( localIndex const p, localIndex const r, localIndex const s )
565  {
566  PsiX[count] = calcModal( p, r, s, X );
567  ++count;
568  } );
569 
570  //Compute the nodal functions at the given points
571  for( int i = 0; i < numNodes; ++i )
572  {
573  N[i] = 0.0;
574  for( int j = 0; j < numNodes; ++j )
575  {
576  N[i] += VDM_inv[i][j] * PsiX[j];
577  }
578  }
579 
580  }
581 
590  static constexpr void calcGradN( localIndex q, arrayView2d< real64 const > VDM_inv, real64 (& gradN)[numNodes][3] )
591  {
592 
593  real64 gradModal[numNodes][3] = {{}};
594  real64 coords[numNodes][3] = {};
595 
596  //Generate the coordinates of the support points
597  generatePointsCoordinates( coords );
598 
599  localIndex count = 0;
600  //Compute the modal functions derivatives at the quadrature point
601  generateIndexes( [&]( localIndex const p, localIndex const r, localIndex const s )
602  {
603  calcGradModal( p, r, s, coords[q], gradModal[count] );
604  ++count;
605  } );
606 
607 
608  //Compute the nodal functions derivatives at the quadrature point
609  for( int i = 0; i < numNodes; ++i )
610  {
611  gradN[i][0] = 0.0;
612  gradN[i][1] = 0.0;
613  gradN[i][2] = 0.0;
614  for( int j = 0; j < numNodes; ++j )
615  {
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];
619  }
620  }
621 
622  }
623 
632  static constexpr void calcGradN( real64 const (&X)[3], arrayView2d< real64 const > VDM_inv, real64 (& gradN)[numNodes][3] )
633  {
634 
635  real64 gradModal[numNodes][3] = {{}};
636 
637  localIndex count = 0;
638  //Compute the modal functions derivatives at the given point
639  generateIndexes( [&]( localIndex const p, localIndex const r, localIndex const s )
640  {
641  calcGradModal( p, r, s, X, gradModal[count] );
642  ++count;
643  } );
644 
645  //Compute the nodal functions derivatives at the given points
646  for( int i = 0; i < numNodes; ++i )
647  {
648  gradN[i][0] = 0.0;
649  gradN[i][1] = 0.0;
650  gradN[i][2] = 0.0;
651  for( int j = 0; j < numNodes; ++j )
652  {
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];
656  }
657  }
658 
659  }
660 
668  template< typename FUNC >
671  static constexpr void computeMassTerm( arrayView2d< real64 const > VDM_inv, FUNC && func )
672  {
673 
674  // Gauss-Legendre points and weights for the quadrature
675  constexpr real64 GLeQuadraturePoints[3] = { -0.7745966692, 0.0, 0.7745966692 };
676  constexpr real64 GLeQuadratureWeights[3] = { 0.5555555556, 0.8888888889, 0.5555555556 };
677  // Gauss-Jacobi points and weights for the quadrature
678  constexpr real64 GJQuadraturePoints[3] = { 0.07299, 0.34700, 0.70500 };
679  constexpr real64 GJQuadratureWeights[3] = { 0.15714, 0.14625, 0.02995 };
680 
681 
682 
683  real64 N[numNodes] = {0.0};
684  for( localIndex a = 0; a < numNodes; ++a )
685  {
686  for( localIndex b = 0; b < numNodes; ++b )
687  {
688  real64 val = 0.0;
689  for( int i = 0; i < 3; ++i )
690  {
691  for( int j = 0; j < 3; ++j )
692  {
693  for( int k = 0; k < 3; ++k )
694  {
695 
696  real64 xi = GLeQuadraturePoints[i];
697  real64 eta = GLeQuadraturePoints[j];
698  real64 chi = GJQuadraturePoints[k];
699 
700  real64 weight = GLeQuadratureWeights[i] * GLeQuadratureWeights[j] * GJQuadratureWeights[k];
701  real64 x_i = (1-chi)*xi;
702  real64 y_j = (1-chi)*eta;
703  real64 z_k = chi;
704  real64 X[3] = {x_i, y_j, z_k};
705 
706  calcN( X, VDM_inv, N );
707 
708  val += weight * N[a] * N[b];
709 
710  }
711  }
712  }
713  func( a, b, val );
714  }
715  }
716  }
717 
725  template< typename FUNC >
728  static void computeStiffnessTerm( arrayView2d< real64 const > VDM_inv, FUNC && func )
729  {
730  // Gauss-Legendre points and weights for the quadrature
731  constexpr real64 GLeQuadraturePoints[3] = { -0.7745966692, 0.0, 0.7745966692 };
732  constexpr real64 GLeQuadratureWeights[3] = { 0.5555555556, 0.8888888889, 0.5555555556 };
733  // Gauss-Jacobi points and weights for the quadrature
734  constexpr real64 GJQuadraturePoints[3] = { 0.07299, 0.34700, 0.70500 };
735  constexpr real64 GJQuadratureWeights[3] = { 0.15714, 0.14625, 0.02995 };
736 
737 
738 
739  real64 gradN[numNodes][3] = {{0.0}};
740 
741  for( localIndex a = 0; a < numNodes; ++a )
742  {
743  for( localIndex b = 0; b < numNodes; ++b )
744  {
745  real64 val = 0.0;
746  for( int i = 0; i < 3; ++i )
747  {
748  for( int j = 0; j < 3; ++j )
749  {
750  for( int k = 0; k < 3; ++k )
751  {
752  real64 xi = GLeQuadraturePoints[i];
753  real64 eta = GLeQuadraturePoints[j];
754  real64 chi = GJQuadraturePoints[k];
755 
756  real64 weight = GLeQuadratureWeights[i] * GLeQuadratureWeights[j] * GJQuadratureWeights[k];
757  real64 x_i = (1-chi)*xi;
758  real64 y_j = (1-chi)*eta;
759  real64 z_k = chi;
760  real64 X[3] = {x_i, y_j, z_k};
761 
762  calcGradN( X, VDM_inv, gradN );
763  real64 dot = gradN[a][0] * gradN[b][0]
764  + gradN[a][1] * gradN[b][1]
765  + gradN[a][2] * gradN[b][2];
766  val += weight * dot;
767  }
768  }
769  }
770  func( a, b, val );
771  }
772  }
773  }
774 
775 };
776 
778 template< int ORDER >
779 class Pk_Pyramid_BCD final : public Pk_Pyramid_BCD_impl< ORDER >, public FiniteElementBase
780 {
781 public:
782 
785 
787  using ImplType::numNodes;
788 
791 
794 
795  Pk_Pyramid_BCD():
799  {}
800 
801 #ifdef __CUDACC__
802  #pragma diag_push
803  #pragma nv_diag_suppress 20012
804 #endif
806  virtual ~Pk_Pyramid_BCD() override final = default;
807 #ifdef __CUDACC__
808  #pragma diag_pop
809 #endif
810 
816  {
817  return static_cast< ImplType * >(this);
818  }
819 
824  ImplType const * getImpl() const
825  {
826  return static_cast< ImplType const * >(this);
827  }
828 
829 };
830 
833 
836 
837 
838 }
839 }
840 
841 #endif // GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_PKPYRAMIDBCD_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:51
Device-compatible (non virtual) Base class for all finite element formulations.
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.
Definition: DataTypes.hpp:191
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:195
Kernel variables (dof numbers, jacobian and residual) located on the stack.