20 #if defined(GEOS_USE_DEVICE) 
   21 #define CALC_FEM_SHAPE_IN_KERNEL 
   26 #ifndef GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_ 
   27 #define GEOS_FINITEELEMENT_ELEMENTFORMULATIONS_FINITEELEMENTBASE_HPP_ 
   31 #include "finiteElement/PDEUtilities.hpp" 
   32 #include "LvArray/src/tensorOps.hpp" 
   62 #ifdef CALC_FEM_SHAPE_IN_KERNEL
 
  109   template< 
typename SUBREGION_TYPE >
 
  122   template< 
typename SUBREGION_TYPE >
 
  126                             SUBREGION_TYPE 
const & cellSubRegion,
 
  146   template< 
typename LEAF, 
typename SUBREGION_TYPE >
 
  150                           SUBREGION_TYPE 
const & cellSubRegion,
 
  154     LEAF::template fillMeshData< SUBREGION_TYPE >( nodeManager,
 
  168   template< 
typename SUBREGION_TYPE >
 
  189   template< 
typename LEAF, 
typename SUBREGION_TYPE >
 
  193               typename LEAF::StackVariables & stack )
 const 
  195     LEAF::setupStack( cellIndex, meshData, stack );
 
  237   template< 
typename LEAF >
 
  241     return LEAF::getNumSupportPoints( stack );
 
  264   template< 
typename LEAF >
 
  268                    real64 const (&X)[LEAF::maxSupportPoints][3],
 
  269                    real64 ( &gradN )[LEAF::maxSupportPoints][3] ) 
const;
 
  283   template< 
typename LEAF >
 
  287                    real64 const (&X)[LEAF::maxSupportPoints][3],
 
  288                    typename LEAF::StackVariables 
const & stack,
 
  289                    real64 ( &gradN )[LEAF::maxSupportPoints][3] ) 
const;
 
  302   template< 
typename LEAF >
 
  307                    real64 ( &gradN )[LEAF::maxSupportPoints][3] ) 
const;
 
  320   template< 
typename LEAF >
 
  325                    typename LEAF::StackVariables 
const & stack,
 
  326                    real64 ( &gradN )[LEAF::maxSupportPoints][3] ) 
const;
 
  339   template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS, 
bool UPPER >
 
  344                                         [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT]
 
  345                                         [MAXSUPPORTPOINTS * NUMDOFSPERTRIALSUPPORTPOINT],
 
  346                                         real64 const & scaleFactor )
 
  363   template< 
typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT, 
bool UPPER = false >
 
  367                                        [LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT]
 
  368                                        [LEAF::maxSupportPoints * NUMDOFSPERTRIALSUPPORTPOINT],
 
  369                                        real64 const scaleFactor = 1.0 )
 const 
  372                                              LEAF::maxSupportPoints,
 
  391   template< localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS >
 
  395                                                  real64 const ( &dofs )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
 
  396                                                  real64 ( & targetVector )[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
 
  397                                                  real64 const scaleFactor )
 
  418   template< 
typename LEAF, localIndex NUMDOFSPERTRIALSUPPORTPOINT >
 
  423                                            real64 const ( &dofs )[LEAF::maxSupportPoints]
 
  424                                            [NUMDOFSPERTRIALSUPPORTPOINT],
 
  425                                            real64 ( & targetVector )[LEAF::maxSupportPoints]
 
  426                                            [NUMDOFSPERTRIALSUPPORTPOINT],
 
  427                                            real64 const scaleFactor = 1.0 )
 const 
  430     addEvaluatedGradGradStabilization< NUMDOFSPERTRIALSUPPORTPOINT, LEAF::maxSupportPoints >( stack,
 
  460   template< 
int NUM_SUPPORT_POINTS >
 
  464               real64 const (&var)[NUM_SUPPORT_POINTS],
 
  472   template< 
int NUM_SUPPORT_POINTS,
 
  477               real64 const (&var)[NUM_SUPPORT_POINTS][NUM_COMPONENTS],
 
  505   template< 
int NUM_SUPPORT_POINTS,
 
  506             typename GRADIENT_TYPE >
 
  509                                  real64 const (&var)[NUM_SUPPORT_POINTS][3],
 
  524   template< 
int NUM_SUPPORT_POINTS,
 
  525             typename GRADIENT_TYPE >
 
  528                                         real64 const (&var)[NUM_SUPPORT_POINTS][3] );
 
  546   template< 
int NUM_SUPPORT_POINTS,
 
  547             typename GRADIENT_TYPE >
 
  550                         real64 const (&var)[NUM_SUPPORT_POINTS],
 
  563   template< 
int NUM_SUPPORT_POINTS,
 
  564             typename GRADIENT_TYPE >
 
  567                         real64 const (&var)[NUM_SUPPORT_POINTS][3],
 
  568                         real64 ( &gradVar )[3][3] );
 
  589   template< 
int NUM_SUPPORT_POINTS,
 
  590             typename GRADIENT_TYPE >
 
  593                                 GRADIENT_TYPE 
const & gradN,
 
  594                                 real64 const (&var)[NUM_SUPPORT_POINTS],
 
  628   template< 
int NUM_SUPPORT_POINTS,
 
  629             typename GRADIENT_TYPE >
 
  632                               real64 const (&var_detJxW)[6],
 
  633                               real64 ( &R )[NUM_SUPPORT_POINTS][3] );
 
  640   template< 
int NUM_SUPPORT_POINTS,
 
  641             typename GRADIENT_TYPE >
 
  644                               real64 const (&var_detJxW)[3][3],
 
  645                               real64 ( &R )[NUM_SUPPORT_POINTS][3] );
 
  655   template< 
int NUM_SUPPORT_POINTS >
 
  658                         real64 const (&forcingTerm_detJxW)[3],
 
  659                         real64 ( &R )[NUM_SUPPORT_POINTS][3] );
 
  678   template< 
int NUM_SUPPORT_POINTS,
 
  679             typename GRADIENT_TYPE >
 
  682                                       real64 const (&var_detJxW)[3][3],
 
  683                                       real64 const (&N)[NUM_SUPPORT_POINTS],
 
  684                                       real64 const (&forcingTerm_detJxW)[3],
 
  685                                       real64 ( &R )[NUM_SUPPORT_POINTS][3] );
 
  692   template< 
int NUM_SUPPORT_POINTS,
 
  693             typename GRADIENT_TYPE >
 
  696                                       real64 const (&var_detJxW)[6],
 
  697                                       real64 const (&N)[NUM_SUPPORT_POINTS],
 
  698                                       real64 const (&forcingTerm_detJxW)[3],
 
  699                                       real64 ( &R )[NUM_SUPPORT_POINTS][3] );
 
  710                           "2nd-dimension of gradN array does not match number of quadrature points" );
 
  713                           "3rd-dimension of gradN array does not match number of support points" );
 
  716                           "4th-dimension of gradN array does not match 3" );
 
  729                           "2nd-dimension of gradN array does not match number of quadrature points" );
 
  773     return PDEUtilities::FunctionSpace::H1;
 
  778 struct FiniteElementBase::FunctionSpaceHelper< 3 >
 
  783     return PDEUtilities::FunctionSpace::H1vector;
 
  791   return FunctionSpaceHelper< N >::getFunctionSpace();
 
  794 template< 
typename LEAF >
 
  799                                     real64 const (&X)[LEAF::maxSupportPoints][3],
 
  800                                     real64 (& gradN)[LEAF::maxSupportPoints][3] )
 const 
  803   return LEAF::calcGradN( q, X, gradN );
 
  806 template< 
typename LEAF >
 
  811                                     real64 const (&X)[LEAF::maxSupportPoints][3],
 
  812                                     typename LEAF::StackVariables 
const & stack,
 
  813                                     real64 ( & gradN )[LEAF::maxSupportPoints][3] )
 const 
  816   return LEAF::calcGradN( q, X, stack, gradN );
 
  819 template< 
typename LEAF >
 
  825                                     real64 (& gradN)[LEAF::maxSupportPoints][3] )
 const 
  829   LvArray::tensorOps::copy< LEAF::maxSupportPoints, 3 >( gradN, 
m_viewGradN[ k ][ q ] );
 
  834 template< 
typename LEAF >
 
  840                                     typename LEAF::StackVariables 
const & stack,
 
  841                                     real64 (& gradN)[LEAF::maxSupportPoints][3] )
 const 
  846   LvArray::tensorOps::copy< LEAF::maxSupportPoints, 3 >( gradN, 
m_viewGradN[ k ][ q ] );
 
  855 template< 
int NUM_SUPPORT_POINTS >
 
  859                                real64 const (&var)[NUM_SUPPORT_POINTS],
 
  862   value = LvArray::tensorOps::AiBi< NUM_SUPPORT_POINTS >( N, var );
 
  865 template< 
int NUM_SUPPORT_POINTS,
 
  870                                real64 const (&var)[NUM_SUPPORT_POINTS][NUM_COMPONENTS],
 
  871                                real64 (& value)[NUM_COMPONENTS] )
 
  874   LvArray::tensorOps::Ri_eq_AjiBj< 3, NUM_SUPPORT_POINTS >( 
value, var, N );
 
  882 template< 
int NUM_SUPPORT_POINTS,
 
  883           typename GRADIENT_TYPE >
 
  887                                            real64 const (&var)[NUM_SUPPORT_POINTS][3],
 
  890   gradVar[0] = gradN[0][0] * var[0][0];
 
  891   gradVar[1] = gradN[0][1] * var[0][1];
 
  892   gradVar[2] = gradN[0][2] * var[0][2];
 
  893   gradVar[3] = gradN[0][2] * var[0][1] + gradN[0][1] * var[0][2];
 
  894   gradVar[4] = gradN[0][2] * var[0][0] + gradN[0][0] * var[0][2];
 
  895   gradVar[5] = gradN[0][1] * var[0][0] + gradN[0][0] * var[0][1];
 
  897   for( 
int a=1; a<NUM_SUPPORT_POINTS; ++a )
 
  899     gradVar[0] = gradVar[0] + gradN[a][0] * var[ a ][0];
 
  900     gradVar[1] = gradVar[1] + gradN[a][1] * var[ a ][1];
 
  901     gradVar[2] = gradVar[2] + gradN[a][2] * var[ a ][2];
 
  902     gradVar[3] = gradVar[3] + gradN[a][2] * var[ a ][1] + gradN[a][1] * var[ a ][2];
 
  903     gradVar[4] = gradVar[4] + gradN[a][2] * var[ a ][0] + gradN[a][0] * var[ a ][2];
 
  904     gradVar[5] = gradVar[5] + gradN[a][1] * var[ a ][0] + gradN[a][0] * var[ a ][1];
 
  908 template< 
int NUM_SUPPORT_POINTS,
 
  909           typename GRADIENT_TYPE >
 
  913                                                   real64 const (&var)[NUM_SUPPORT_POINTS][3] )
 
  915   real64 result = gradN[0][0] * var[0][0] + gradN[0][1] * var[0][1] + gradN[0][2] * var[0][2];
 
  917   for( 
int a=1; a<NUM_SUPPORT_POINTS; ++a )
 
  919     result = result + gradN[a][0] * var[a][0] + gradN[a][1] * var[a][1] + gradN[a][2] * var[a][2];
 
  924 template< 
int NUM_SUPPORT_POINTS,
 
  925           typename GRADIENT_TYPE >
 
  929                                   real64 const (&var)[NUM_SUPPORT_POINTS],
 
  932   LvArray::tensorOps::Ri_eq_AjiBj< 3, NUM_SUPPORT_POINTS >( gradVar, gradN, var );
 
  935 template< 
int NUM_SUPPORT_POINTS,
 
  936           typename GRADIENT_TYPE >
 
  940                                   real64 const (&var)[NUM_SUPPORT_POINTS][3],
 
  941                                   real64 (& gradVar)[3][3] )
 
  943   LvArray::tensorOps::Rij_eq_AkiBkj< 3, 3, NUM_SUPPORT_POINTS >( gradVar, var, gradN );
 
  948 template< 
int NUM_SUPPORT_POINTS,
 
  949           typename GRADIENT_TYPE >
 
  953                                           GRADIENT_TYPE 
const & gradN,
 
  954                                           real64 const (&var)[NUM_SUPPORT_POINTS],
 
  958   value = N[0] * var[0];
 
  959   for( 
int i = 0; i < 3; ++i )
 
  961     gradVar[i] = var[0] * gradN[0][i];
 
  964   for( 
int a=1; a<NUM_SUPPORT_POINTS; ++a )
 
  967     for( 
int i = 0; i < 3; ++i )
 
  969       gradVar[i] = gradVar[i] + var[ a ] * gradN[a][i];
 
  976 template< 
int NUM_SUPPORT_POINTS,
 
  977           typename GRADIENT_TYPE >
 
  981                                         real64 const (&var_detJxW)[6],
 
  982                                         real64 (& R)[NUM_SUPPORT_POINTS][3] )
 
  984   for( 
int a=0; a<NUM_SUPPORT_POINTS; ++a )
 
  986     R[a][0] = R[a][0] + var_detJxW[0] * gradN[a][0] + var_detJxW[5] * gradN[a][1] + var_detJxW[4] * gradN[a][2];
 
  987     R[a][1] = R[a][1] + var_detJxW[5] * gradN[a][0] + var_detJxW[1] * gradN[a][1] + var_detJxW[3] * gradN[a][2];
 
  988     R[a][2] = R[a][2] + var_detJxW[4] * gradN[a][0] + var_detJxW[3] * gradN[a][1] + var_detJxW[2] * gradN[a][2];
 
  993 template< 
int NUM_SUPPORT_POINTS,
 
  994           typename GRADIENT_TYPE >
 
  998                                         real64 const (&var_detJxW)[3][3],
 
  999                                         real64 (& R)[NUM_SUPPORT_POINTS][3] )
 
 1001   for( 
int a=0; a<NUM_SUPPORT_POINTS; ++a )
 
 1003     LvArray::tensorOps::Ri_add_AijBj< 3, 3 >( R[a], var_detJxW, gradN[a] );
 
 1007 template< 
int NUM_SUPPORT_POINTS >
 
 1011                                   real64 const (&var_detJxW)[3],
 
 1012                                   real64 ( & R )[NUM_SUPPORT_POINTS][3] )
 
 1014   for( 
int a=0; a<NUM_SUPPORT_POINTS; ++a )
 
 1016     LvArray::tensorOps::scaledAdd< 3 >( R[a], var_detJxW, N[a] );
 
 1021 template< 
int NUM_SUPPORT_POINTS,
 
 1022           typename GRADIENT_TYPE >
 
 1026                                                 real64 const (&var_detJxW)[6],
 
 1027                                                 real64 const (&N)[NUM_SUPPORT_POINTS],
 
 1028                                                 real64 const (&forcingTerm_detJxW)[3],
 
 1029                                                 real64 (& R)[NUM_SUPPORT_POINTS][3] )
 
 1031   for( 
int a=0; a<NUM_SUPPORT_POINTS; ++a )
 
 1033     R[a][0] = R[a][0] + var_detJxW[0] * gradN[a][0] + var_detJxW[5] * gradN[a][1] + var_detJxW[4] * gradN[a][2] + forcingTerm_detJxW[0] * N[a];
 
 1034     R[a][1] = R[a][1] + var_detJxW[5] * gradN[a][0] + var_detJxW[1] * gradN[a][1] + var_detJxW[3] * gradN[a][2] + forcingTerm_detJxW[1] * N[a];
 
 1035     R[a][2] = R[a][2] + var_detJxW[4] * gradN[a][0] + var_detJxW[3] * gradN[a][1] + var_detJxW[2] * gradN[a][2] + forcingTerm_detJxW[2] * N[a];
 
 1039 template< 
int NUM_SUPPORT_POINTS,
 
 1040           typename GRADIENT_TYPE >
 
 1044                                                 real64 const (&var_detJxW)[3][3],
 
 1045                                                 real64 const (&N)[NUM_SUPPORT_POINTS],
 
 1046                                                 real64 const (&forcingTerm_detJxW)[3],
 
 1047                                                 real64 (& R)[NUM_SUPPORT_POINTS][3] )
 
 1049   for( 
int a=0; a<NUM_SUPPORT_POINTS; ++a )
 
 1051     R[a][0] = R[a][0] + var_detJxW[0][0] * gradN[a][0] + var_detJxW[0][1] * gradN[a][1] + var_detJxW[0][2] * gradN[a][2] + forcingTerm_detJxW[0] * N[a];
 
 1052     R[a][1] = R[a][1] + var_detJxW[1][0] * gradN[a][0] + var_detJxW[1][1] * gradN[a][1] + var_detJxW[1][2] * gradN[a][2] + forcingTerm_detJxW[1] * N[a];
 
 1053     R[a][2] = R[a][2] + var_detJxW[2][0] * gradN[a][0] + var_detJxW[2][1] * gradN[a][1] + var_detJxW[2][2] * gradN[a][2] + forcingTerm_detJxW[2] * N[a];
 
 1063 #define USING_FINITEELEMENTBASE                       \ 
 1064   using FiniteElementBase::value;                     \ 
 1065   using FiniteElementBase::symmetricGradient;         \ 
 1066   using FiniteElementBase::gradient;                  \ 
 1067   using FiniteElementBase::valueAndGradient;          \ 
 1068   using FiniteElementBase::plusGradNajAij;           \ 
 1069   using FiniteElementBase::plusNaFi;                 \ 
 1070   using FiniteElementBase::plusGradNajAijPlusNaFi; 
#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_ERROR_IF_NE_MSG(lhs, rhs, msg)
Raise a hard error if two values are not equal.
 
This class provides an interface to ObjectManagerBase in order to manage edge data.
 
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
 
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
 
Base class for FEM element implementations.
 
static GEOS_HOST_DEVICE void value(real64 const (&N)[NUM_SUPPORT_POINTS], real64 const (&var)[NUM_SUPPORT_POINTS][NUM_COMPONENTS], real64(&value)[NUM_COMPONENTS])
Compute the interpolated value of a vector variable.
 
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints() const =0
Virtual getter for the number of support points per element.
 
static GEOS_HOST_DEVICE void plusGradNajAijPlusNaFi(GRADIENT_TYPE const &gradN, real64 const (&var_detJxW)[3][3], real64 const (&N)[NUM_SUPPORT_POINTS], real64 const (&forcingTerm_detJxW)[3], real64(&R)[NUM_SUPPORT_POINTS][3])
Inner product of each basis function gradient with a rank-2 symmetric tensor added to the product eac...
 
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints() const =0
Get the maximum number of support points for this element.
 
static void fillMeshData(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, MeshData< SUBREGION_TYPE > &meshData)
Method to fill a MeshData object.
 
arrayView2d< real64 const > getDetJView() const
Getter for m_viewDetJ.
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void setupStack(localIndex const &cellIndex, MeshData< SUBREGION_TYPE > const &meshData, StackVariables &stack)
Empty setup method.
 
static GEOS_HOST_DEVICE void plusGradNajAijPlusNaFi(GRADIENT_TYPE const &gradN, real64 const (&var_detJxW)[6], real64 const (&N)[NUM_SUPPORT_POINTS], real64 const (&forcingTerm_detJxW)[3], real64(&R)[NUM_SUPPORT_POINTS][3])
Inner product of each basis function gradient with a rank-2 tensor added to the product each shape fu...
 
virtual GEOS_HOST_DEVICE ~FiniteElementBase()
Destructor.
 
static GEOS_HOST_DEVICE void gradient(GRADIENT_TYPE const &gradN, real64 const (&var)[NUM_SUPPORT_POINTS], real64(&gradVar)[3])
Calculate the gradient of a scalar valued support field at a point using the input basis function gra...
 
static GEOS_HOST_DEVICE void plusGradNajAij(GRADIENT_TYPE const &gradN, real64 const (&var_detJxW)[3][3], real64(&R)[NUM_SUPPORT_POINTS][3])
Inner product of each basis function gradient with a rank-2 symmetric tensor.
 
FiniteElementBase & operator=(FiniteElementBase const &)=delete
Deleted copy assignment operator.
 
GEOS_HOST_DEVICE void setup(localIndex const &cellIndex, typename LEAF::template MeshData< SUBREGION_TYPE > const &meshData, typename LEAF::StackVariables &stack) const
Abstract setup method, possibly computing cell-dependent properties.
 
arrayView4d< real64 const > m_viewGradN
View to potentially hold pre-calculated shape function gradients.
 
GEOS_HOST_DEVICE FiniteElementBase(FiniteElementBase const &source)
Copy Constructor.
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addGradGradStabilization(StackVariables const &stack, real64(&matrix)[MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT], real64 const &scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
 
GEOS_HOST_DEVICE real64 getGradN(localIndex const k, localIndex const q, real64 const (&X)[LEAF::maxSupportPoints][3], typename LEAF::StackVariables const &stack, real64(&gradN)[LEAF::maxSupportPoints][3]) const
Get the shape function gradients.
 
void setDetJView(arrayView2d< real64 const > const &source)
Sets m_viewDetJ equal to an input view.
 
static GEOS_HOST_DEVICE void plusGradNajAij(GRADIENT_TYPE const &gradN, real64 const (&var_detJxW)[6], real64(&R)[NUM_SUPPORT_POINTS][3])
Inner product of each basis function gradient with a rank-2 symmetric tensor.
 
FiniteElementBase()=default
Default Constructor.
 
GEOS_HOST_DEVICE real64 getGradN(localIndex const k, localIndex const q, int const X, real64(&gradN)[LEAF::maxSupportPoints][3]) const
Get the shape function gradients.
 
GEOS_HOST_DEVICE localIndex numSupportPoints(typename LEAF::StackVariables const &stack) const
Getter for the number of support points per element.
 
static GEOS_HOST_DEVICE void symmetricGradient(GRADIENT_TYPE const &gradN, real64 const (&var)[NUM_SUPPORT_POINTS][3], real64(&gradVar)[6])
Calculate the symmetric gradient of a vector valued support field at a point using the stored basis f...
 
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints() const =0
Virtual getter for the number of quadrature points per element.
 
static void initialize(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &cellSubRegion, typename LEAF::template MeshData< SUBREGION_TYPE > &meshData)
Abstract initialization method.
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void addEvaluatedGradGradStabilization(StackVariables const &stack, real64 const (&dofs)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor)
Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilin...
 
FiniteElementBase & operator=(FiniteElementBase &&)=delete
Deleted move assignment operator.
 
static GEOS_HOST_DEVICE real64 symmetricGradientTrace(GRADIENT_TYPE const &gradN, real64 const (&var)[NUM_SUPPORT_POINTS][3])
Calculate the trace of the symmetric gradient of a vector valued support field (i....
 
GEOS_HOST_DEVICE void addGradGradStabilizationMatrix(typename LEAF::StackVariables const &stack, real64(&matrix)[LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add stabilization of grad-grad bilinear form to input matrix.
 
static GEOS_HOST_DEVICE void plusNaFi(real64 const (&N)[NUM_SUPPORT_POINTS], real64 const (&forcingTerm_detJxW)[3], real64(&R)[NUM_SUPPORT_POINTS][3])
Product of each shape function with a vector forcing term.
 
constexpr static GEOS_HOST_DEVICE PDEUtilities::FunctionSpace getFunctionSpace()
Getter for the function space.
 
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void addEvaluatedGradGradStabilizationVector(typename LEAF::StackVariables const &stack, real64 const (&dofs)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64(&targetVector)[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT], real64 const scaleFactor=1.0) const
Add a grad-grad stabilization operator evaluated at a provided vector of dofs to input vector.
 
GEOS_HOST_DEVICE real64 getGradN(localIndex const k, localIndex const q, real64 const (&X)[LEAF::maxSupportPoints][3], real64(&gradN)[LEAF::maxSupportPoints][3]) const
Get the shape function gradients.
 
static GEOS_HOST_DEVICE void value(real64 const (&N)[NUM_SUPPORT_POINTS], real64 const (&var)[NUM_SUPPORT_POINTS], real64 &value)
Compute the interpolated value of a variable.
 
constexpr static int numSamplingPointsPerDirection
Number of sampling points.
 
arrayView4d< real64 const > getGradNView() const
Getter for m_viewGradN.
 
FiniteElementBase(FiniteElementBase &&)=default
Default Move constructor.
 
static GEOS_HOST_DEVICE void gradient(GRADIENT_TYPE const &gradN, real64 const (&var)[NUM_SUPPORT_POINTS][3], real64(&gradVar)[3][3])
Calculate the gradient of a vector valued support field at a point using the input basis function gra...
 
GEOS_HOST_DEVICE real64 getGradN(localIndex const k, localIndex const q, int const X, typename LEAF::StackVariables const &stack, real64(&gradN)[LEAF::maxSupportPoints][3]) const
Get the shape function gradients.
 
arrayView2d< real64 const > m_viewDetJ
 
void setGradNView(arrayView4d< real64 const > const &source)
Sets m_viewGradN equal to an input view.
 
static GEOS_HOST_DEVICE void valueAndGradient(real64 const (&N)[NUM_SUPPORT_POINTS], GRADIENT_TYPE const &gradN, real64 const (&var)[NUM_SUPPORT_POINTS], real64 &value, real64(&gradVar)[3])
Calculate the value and gradient of a scalar valued support field at a point using the input basis fu...
 
double real64
64-bit floating point type.
 
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
 
ArrayView< T, 4, USD > arrayView4d
Alias for 4D array view.
 
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
 
An helper struct to determine the function space.
 
Variables used to initialize the class.
 
Kernel variables allocated on the stack.