GEOS
Classes | Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
geos::finiteElement::FiniteElementBase Class Referenceabstract

Base class for FEM element implementations. More...

#include <FiniteElementBase.hpp>

Inheritance diagram for geos::finiteElement::FiniteElementBase:
Inheritance graph
[legend]

Classes

struct  FunctionSpaceHelper
 An helper struct to determine the function space. More...
 
struct  MeshData
 Variables used to initialize the class. More...
 
struct  StackVariables
 Kernel variables allocated on the stack. More...
 

Public Member Functions

 FiniteElementBase ()=default
 Default Constructor.
 
GEOS_HOST_DEVICE FiniteElementBase (FiniteElementBase const &source)
 Copy Constructor. More...
 
 FiniteElementBase (FiniteElementBase &&)=default
 Default Move constructor.
 
FiniteElementBaseoperator= (FiniteElementBase const &)=delete
 Deleted copy assignment operator. More...
 
FiniteElementBaseoperator= (FiniteElementBase &&)=delete
 Deleted move assignment operator. More...
 
virtual GEOS_HOST_DEVICE ~FiniteElementBase ()
 Destructor.
 
template<typename LEAF , typename SUBREGION_TYPE >
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. More...
 
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints () const =0
 Virtual getter for the number of quadrature points per element. More...
 
virtual GEOS_HOST_DEVICE localIndex getNumSupportPoints () const =0
 Virtual getter for the number of support points per element. More...
 
template<typename LEAF >
GEOS_HOST_DEVICE localIndex numSupportPoints (typename LEAF::StackVariables const &stack) const
 Getter for the number of support points per element. More...
 
virtual GEOS_HOST_DEVICE localIndex getMaxSupportPoints () const =0
 Get the maximum number of support points for this element. More...
 
template<typename LEAF >
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. More...
 
template<typename LEAF >
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. More...
 
template<typename LEAF >
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. More...
 
template<typename LEAF >
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. More...
 
template<typename LEAF , localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER = false>
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. More...
 
template<typename LEAF , localIndex NUMDOFSPERTRIALSUPPORTPOINT>
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. More...
 

Static Public Member Functions

template<typename SUBREGION_TYPE >
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. More...
 
template<typename LEAF , typename SUBREGION_TYPE >
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. More...
 
template<typename SUBREGION_TYPE >
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void setupStack (localIndex const &cellIndex, MeshData< SUBREGION_TYPE > const &meshData, StackVariables &stack)
 Empty setup method. More...
 
template<int N>
constexpr static GEOS_HOST_DEVICE PDEUtilities::FunctionSpace getFunctionSpace ()
 Getter for the function space. More...
 
template<localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS, bool UPPER>
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 bilinear form. More...
 
template<localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS>
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 bilinear form. More...
 
Value Operator Functions
template<int NUM_SUPPORT_POINTS>
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. More...
 
template<int NUM_SUPPORT_POINTS, int NUM_COMPONENTS>
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. More...
 
Gradient Operator Functions
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 function gradients for all support points. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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.e. the volumetric strain for the displacement field) at a point using the stored basis function gradients for all support points. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 gradients. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 gradients. More...
 
Multi-Operator Functions
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 function gradients. More...
 

Static Public Attributes

constexpr static int numSamplingPointsPerDirection = 10
 Number of sampling points.
 

Scattering Operator Functions

These functions take quadrature data and map it to the support points through some operator.

arrayView4d< real64 const > m_viewGradN
 View to potentially hold pre-calculated shape function gradients.
 
arrayView2d< real64 const > m_viewDetJ
 
void setGradNView (arrayView4d< real64 const > const &source)
 Sets m_viewGradN equal to an input view. More...
 
void setDetJView (arrayView2d< real64 const > const &source)
 Sets m_viewDetJ equal to an input view. More...
 
arrayView4d< real64 const > getGradNView () const
 Getter for m_viewGradN. More...
 
arrayView2d< real64 const > getDetJView () const
 Getter for m_viewDetJ. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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. More...
 
template<int NUM_SUPPORT_POINTS>
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. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 each shape function with a vector. More...
 
template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
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 function with a vector. More...
 

Detailed Description

Base class for FEM element implementations.

Definition at line 45 of file FiniteElementBase.hpp.

Constructor & Destructor Documentation

◆ FiniteElementBase()

GEOS_HOST_DEVICE geos::finiteElement::FiniteElementBase::FiniteElementBase ( FiniteElementBase const &  source)
inline

Copy Constructor.

Parameters
sourceThe object to copy.

Definition at line 61 of file FiniteElementBase.hpp.

Member Function Documentation

◆ addEvaluatedGradGradStabilization()

template<localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::FiniteElementBase::addEvaluatedGradGradStabilization ( StackVariables const &  stack,
real64 const (&)  dofs[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
real64(&)  targetVector[MAXSUPPORTPOINTS][NUMDOFSPERTRIALSUPPORTPOINT],
real64 const  scaleFactor 
)
inlinestatic

Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilinear form.

This method is intended to be used with targetVector being the residual and dofs being the degrees of freedom of the previous solution.

Template Parameters
NUMDOFSPERTRIALSUPPORTPOINTNumber of degrees of freedom for each support point.
Parameters
stackStack variables as filled by setupStack.
dofsThe degrees of freedom of the function where the stabilization operator has to be evaluated.
targetVectorThe input vector to which values have to be added, seen in chunks of length NUMDOFSPERTRIALSUPPORTPOINT.
scaleFactorScaling of the stabilization matrix.

Definition at line 407 of file FiniteElementBase.hpp.

◆ addEvaluatedGradGradStabilizationVector()

template<typename LEAF , localIndex NUMDOFSPERTRIALSUPPORTPOINT>
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void geos::finiteElement::FiniteElementBase::addEvaluatedGradGradStabilizationVector ( typename LEAF::StackVariables const &  stack,
real64 const (&)  dofs[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT],
real64(&)  targetVector[LEAF::maxSupportPoints][NUMDOFSPERTRIALSUPPORTPOINT],
real64 const  scaleFactor = 1.0 
) const
inline

Add a grad-grad stabilization operator evaluated at a provided vector of dofs to input vector.

This method is used to modify a residual consistently when the jacobian includes a stabilization term.

Template Parameters
LEAFType of the derived finite element implementation.
NUMDOFSPERTRIALSUPPORTPOINTNumber of degrees of freedom for each support point.
Parameters
stackStack variables created by a call to setup.
dofsThe vector of dofs to evaluate the stabilization.
targetVectorThe input vector to which values have to be added, seen in chunks of length NUMDOFSPERTRIALSUPPORTPOINT.
scaleFactorOptional scaling of the stabilization matrix. Defaults to 1.0.

Definition at line 435 of file FiniteElementBase.hpp.

◆ addGradGradStabilization()

template<localIndex NUMDOFSPERTRIALSUPPORTPOINT, localIndex MAXSUPPORTPOINTS, bool UPPER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::FiniteElementBase::addGradGradStabilization ( StackVariables const &  stack,
real64(&)  matrix[MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT][MAXSUPPORTPOINTS *NUMDOFSPERTRIALSUPPORTPOINT],
real64 const &  scaleFactor 
)
inlinestatic

Empty method, here for compatibility with methods that require a stabilization of the grad-grad bilinear form.

Template Parameters
NUMDOFSPERTRIALSUPPORTPOINTNumber of degrees of freedom for each support point.
MAXSUPPORTPOINTSMaximum number of support points allowed for this element.
UPPERIf true only the upper triangular part of matrix is modified.
Parameters
stackStack variables as filled by setupStack.
matrixThe matrix that needs to be stabilized.
scaleFactorScaling of the stabilization matrix.

Definition at line 355 of file FiniteElementBase.hpp.

◆ addGradGradStabilizationMatrix()

template<typename LEAF , localIndex NUMDOFSPERTRIALSUPPORTPOINT, bool UPPER = false>
GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::addGradGradStabilizationMatrix ( typename LEAF::StackVariables const &  stack,
real64(&)  matrix[LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT][LEAF::maxSupportPoints *NUMDOFSPERTRIALSUPPORTPOINT],
real64 const  scaleFactor = 1.0 
) const
inline

Add stabilization of grad-grad bilinear form to input matrix.

Template Parameters
LEAFType of the derived finite element implementation.
NUMDOFSPERTRIALSUPPORTPOINTNumber of degrees of freedom for each support point.
UPPERIf true only the upper triangular part of matrix is modified.
Parameters
stackStack variables created by a call to setup.
matrixThe input matrix to which values have to be added.
scaleFactorOptional scaling of the stabilization matrix. Defaults to 1.0.

Definition at line 378 of file FiniteElementBase.hpp.

◆ fillMeshData()

template<typename SUBREGION_TYPE >
static void geos::finiteElement::FiniteElementBase::fillMeshData ( NodeManager const &  nodeManager,
EdgeManager const &  edgeManager,
FaceManager const &  faceManager,
SUBREGION_TYPE const &  cellSubRegion,
MeshData< SUBREGION_TYPE > &  meshData 
)
inlinestatic

Method to fill a MeshData object.

Parameters
nodeManagerThe node manager.
edgeManagerThe edge manager.
faceManagerThe face manager.
cellSubRegionThe cell sub-region for which the element has to be initialized.
meshDataMeshData struct to be filled.

Definition at line 136 of file FiniteElementBase.hpp.

◆ getDetJView()

arrayView2d< real64 const > geos::finiteElement::FiniteElementBase::getDetJView ( ) const
inline

Getter for m_viewDetJ.

Returns
A new arrayView copy of m_viewDetJ.

Definition at line 759 of file FiniteElementBase.hpp.

◆ getFunctionSpace()

template<int N>
constexpr static GEOS_HOST_DEVICE PDEUtilities::FunctionSpace geos::finiteElement::FiniteElementBase::getFunctionSpace ( )
staticconstexpr

Getter for the function space.

Template Parameters
Thenumber of components per support point (i.e., 1 if scalar variable, 3 if vector variable)
Returns
The function space.

◆ getGradN() [1/4]

template<typename LEAF >
GEOS_HOST_DEVICE real64 geos::finiteElement::FiniteElementBase::getGradN ( localIndex const  k,
localIndex const  q,
int const  X,
real64(&)  gradN[LEAF::maxSupportPoints][3] 
) const

Get the shape function gradients.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
kThe element index.
qThe quadrature point index.
Xdummy variable.
gradNReturn array of the shape function gradients.
Returns
The determinant of the Jacobian transformation matrix.

This function returns pre-calculated shape function gradients.

◆ getGradN() [2/4]

template<typename LEAF >
GEOS_HOST_DEVICE real64 geos::finiteElement::FiniteElementBase::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.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
kThe element index.
qThe quadrature point index.
Xdummy variable.
stackStack variables relative to the element k created by a call to setup.
gradNReturn array of the shape function gradients.
Returns
The determinant of the Jacobian transformation matrix.

This function returns pre-calculated shape function gradients.

◆ getGradN() [3/4]

template<typename LEAF >
GEOS_HOST_DEVICE real64 geos::finiteElement::FiniteElementBase::getGradN ( localIndex const  k,
localIndex const  q,
real64 const (&)  X[LEAF::maxSupportPoints][3],
real64(&)  gradN[LEAF::maxSupportPoints][3] 
) const

Get the shape function gradients.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
kThe element index.
qThe quadrature point index.
XArray of coordinates as the reference for the gradients.
gradNReturn array of the shape function gradients.
Returns
The determinant of the Jacobian transformation matrix.

This function calls the function to calculate shape function gradients.

◆ getGradN() [4/4]

template<typename LEAF >
GEOS_HOST_DEVICE real64 geos::finiteElement::FiniteElementBase::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.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
kThe element index.
qThe quadrature point index.
XArray of coordinates as the reference for the gradients.
stackStack variables relative to the element k created by a call to setup.
gradNReturn array of the shape function gradients.
Returns
The determinant of the Jacobian transformation matrix.

This function calls the function to calculate shape function gradients.

◆ getGradNView()

arrayView4d< real64 const > geos::finiteElement::FiniteElementBase::getGradNView ( ) const
inline

Getter for m_viewGradN.

Returns
A new arrayView copy of m_viewGradN.

Definition at line 750 of file FiniteElementBase.hpp.

◆ getMaxSupportPoints()

virtual GEOS_HOST_DEVICE localIndex geos::finiteElement::FiniteElementBase::getMaxSupportPoints ( ) const
pure virtual

◆ getNumQuadraturePoints()

virtual GEOS_HOST_DEVICE localIndex geos::finiteElement::FiniteElementBase::getNumQuadraturePoints ( ) const
pure virtual

◆ getNumSupportPoints()

virtual GEOS_HOST_DEVICE localIndex geos::finiteElement::FiniteElementBase::getNumSupportPoints ( ) const
pure virtual

◆ gradient() [1/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::gradient ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var[NUM_SUPPORT_POINTS],
real64(&)  gradVar[3] 
)
static

Calculate the gradient of a scalar valued support field at a point using the input basis function gradients.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
varThe vector valued support field that the gradient operator will be applied to.
gradVarThe gradient.

More precisely, the operator is defined as:

\[ grad_{j} = \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{a}\right ), \]

◆ gradient() [2/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::gradient ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var[NUM_SUPPORT_POINTS][3],
real64(&)  gradVar[3][3] 
)
static

Calculate the gradient of a vector valued support field at a point using the input basis function gradients.

Calculate the gradient of a scalar valued support field at a point using the input basis function gradients.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
varThe vector valued support field that the gradient operator will be applied to.
gradVarThe gradient.

More precisely, the operator is defined as:

\[ grad_{j} = \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{a}\right ), \]

More precisely, the operator is defined as:

\[ grad_{ij} = \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{ai}\right ), \]

◆ initialize()

template<typename LEAF , typename SUBREGION_TYPE >
static void geos::finiteElement::FiniteElementBase::initialize ( NodeManager const &  nodeManager,
EdgeManager const &  edgeManager,
FaceManager const &  faceManager,
SUBREGION_TYPE const &  cellSubRegion,
typename LEAF::template MeshData< SUBREGION_TYPE > &  meshData 
)
inlinestatic

Abstract initialization method.

It calls the fillMeshData method of the specific element implementation.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
nodeManagerThe node manager.
edgeManagerThe edge manager.
faceManagerThe face manager.
cellSubRegionThe cell sub-region for which the element has to be initialized.
meshDataThe struct to be filled according to the LEAF class needs.

Definition at line 160 of file FiniteElementBase.hpp.

◆ numSupportPoints()

template<typename LEAF >
GEOS_HOST_DEVICE localIndex geos::finiteElement::FiniteElementBase::numSupportPoints ( typename LEAF::StackVariables const &  stack) const
inline

Getter for the number of support points per element.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
stackStack variables created by a call to setup.
Returns
The number of support points per element.

Definition at line 252 of file FiniteElementBase.hpp.

◆ operator=() [1/2]

FiniteElementBase& geos::finiteElement::FiniteElementBase::operator= ( FiniteElementBase &&  )
delete

Deleted move assignment operator.

Returns
deleted

◆ operator=() [2/2]

FiniteElementBase& geos::finiteElement::FiniteElementBase::operator= ( FiniteElementBase const &  )
delete

Deleted copy assignment operator.

Returns
deleted

◆ plusGradNajAij() [1/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::plusGradNajAij ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var_detJxW[3][3],
real64(&)  R[NUM_SUPPORT_POINTS][3] 
)
static

Inner product of each basis function gradient with a rank-2 symmetric tensor.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
var_detJxWThe rank-2 tensor at q scaled by J*W.
RThe vector at each support point which will hold the result from the tensor contraction.

More precisely, the operator is defined as:

\[ R_i = \sum_a^{nSupport} \left( \frac{\partial N_a}{\partial X_j} var_{ij}\right), \]

where $\frac{\partial N_a}{\partial X_j}$ is the basis function gradient, $var_{ij}$ is the rank-2 symmetric tensor.

Inner product of each basis function gradient with a rank-2 tensor.

◆ plusGradNajAij() [2/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::plusGradNajAij ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var_detJxW[6],
real64(&)  R[NUM_SUPPORT_POINTS][3] 
)
static

Inner product of each basis function gradient with a rank-2 symmetric tensor.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
var_detJxWThe rank-2 tensor at q scaled by J*W.
RThe vector at each support point which will hold the result from the tensor contraction.

More precisely, the operator is defined as:

\[ R_i = \sum_a^{nSupport} \left( \frac{\partial N_a}{\partial X_j} var_{ij}\right), \]

where $\frac{\partial N_a}{\partial X_j}$ is the basis function gradient, $var_{ij}$ is the rank-2 symmetric tensor.

◆ plusGradNajAijPlusNaFi() [1/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::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] 
)
static

Inner product of each basis function gradient with a rank-2 symmetric tensor added to the product each shape function with a vector.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
var_detJxWThe rank-2 symmetric tensor at q scaled by J*W.
NThe shape function value at a predetermined coordinate in the element.
forcingTerm_detJxWA vector scaled by detJxW
RThe vector at each support point which will hold the result from the tensor contraction.

\[ R_i = \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{ij} + N_a f_i \right ), \]

◆ plusGradNajAijPlusNaFi() [2/2]

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::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] 
)
static

Inner product of each basis function gradient with a rank-2 tensor added to the product each shape function with a vector.

Inner product of each basis function gradient with a rank-2 symmetric tensor added to the product each shape function with a vector.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape function gradients.
Parameters
gradNThe basis function gradients at a point in the element.
var_detJxWThe rank-2 symmetric tensor at q scaled by J*W.
NThe shape function value at a predetermined coordinate in the element.
forcingTerm_detJxWA vector scaled by detJxW
RThe vector at each support point which will hold the result from the tensor contraction.

\[ R_i = \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{ij} + N_a f_i \right ), \]

◆ plusNaFi()

template<int NUM_SUPPORT_POINTS>
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::plusNaFi ( real64 const (&)  N[NUM_SUPPORT_POINTS],
real64 const (&)  forcingTerm_detJxW[3],
real64(&)  R[NUM_SUPPORT_POINTS][3] 
)
static

Product of each shape function with a vector forcing term.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
Parameters
NThe shape function value at a predetermined coordinate in the element.
forcingTerm_detJxWA vector scaled by detJxW
RThe vector at each support point which will hold the result from the tensor contraction.

◆ setDetJView()

void geos::finiteElement::FiniteElementBase::setDetJView ( arrayView2d< real64 const > const &  source)
inline

Sets m_viewDetJ equal to an input view.

Parameters
sourceThe view to assign to m_viewDetJ.

Definition at line 738 of file FiniteElementBase.hpp.

◆ setGradNView()

void geos::finiteElement::FiniteElementBase::setGradNView ( arrayView4d< real64 const > const &  source)
inline

Sets m_viewGradN equal to an input view.

Parameters
sourceThe view to assign to m_viewGradN.

Definition at line 719 of file FiniteElementBase.hpp.

◆ setup()

template<typename LEAF , typename SUBREGION_TYPE >
GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::setup ( localIndex const &  cellIndex,
typename LEAF::template MeshData< SUBREGION_TYPE > const &  meshData,
typename LEAF::StackVariables &  stack 
) const
inline

Abstract setup method, possibly computing cell-dependent properties.

Template Parameters
LEAFType of the derived finite element implementation.
Parameters
cellIndexThe index of the cell with respect to the cell sub region to which the element has been initialized previously (see initialize).
meshDataA MeshData object previously filled.
stackObject that holds stack variables.

Definition at line 204 of file FiniteElementBase.hpp.

◆ setupStack()

template<typename SUBREGION_TYPE >
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::FiniteElementBase::setupStack ( localIndex const &  cellIndex,
MeshData< SUBREGION_TYPE > const &  meshData,
StackVariables stack 
)
inlinestatic

Empty setup method.

Parameters
cellIndexThe index of the cell with respect to the cell sub region.
meshDataMeshData struct filled by fillMeshData.
stackObject that holds stack variables.

Definition at line 184 of file FiniteElementBase.hpp.

◆ symmetricGradient()

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::symmetricGradient ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var[NUM_SUPPORT_POINTS][3],
real64(&)  gradVar[6] 
)
static

Calculate the symmetric gradient of a vector valued support field at a point using the stored basis function gradients for all support points.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape
Parameters
gradNThe basis function gradients at a point in the element.
varThe vector valued support field that the gradient operator will be applied to.
gradVarThe symmetric gradient in Voigt notation.

More precisely, the operator is defined as:

\[ grad^s_{ij} = \frac{1}{2} \sum_a^{nSupport} \left ( \frac{\partial N_a}{\partial X_j} var_{ai} + \frac{\partial N_a}{\partial X_i} var_{aj}\right ), \]

◆ symmetricGradientTrace()

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE real64 geos::finiteElement::FiniteElementBase::symmetricGradientTrace ( GRADIENT_TYPE const &  gradN,
real64 const (&)  var[NUM_SUPPORT_POINTS][3] 
)
static

Calculate the trace of the symmetric gradient of a vector valued support field (i.e. the volumetric strain for the displacement field) at a point using the stored basis function gradients for all support points.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape
Parameters
gradNThe basis function gradients at a point in the element.
varThe vector valued support field that the gradient operator will be applied to.
Returns
The trace of the symetric gradient tensor.

◆ value() [1/2]

template<int NUM_SUPPORT_POINTS>
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::value ( real64 const (&)  N[NUM_SUPPORT_POINTS],
real64 const (&)  var[NUM_SUPPORT_POINTS],
real64 value 
)
static

Compute the interpolated value of a variable.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
Parameters
NArray (for each support point) of shape function values at the coordinate the variable is to be interpolated.
varArray of variable values for each support point.
valueThe interpolated value of var.

This is the standard finite element interpolation operator of a discrete variable defined at the support points. The operator is expressed as:

\[ value = \sum_a^{numSupport} \left ( N_a var_a \right ), \]

Note
The shape function values N must be evaluated prior to calling this function.

◆ value() [2/2]

template<int NUM_SUPPORT_POINTS, int NUM_COMPONENTS>
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::value ( real64 const (&)  N[NUM_SUPPORT_POINTS],
real64 const (&)  var[NUM_SUPPORT_POINTS][NUM_COMPONENTS],
real64(&)  value[NUM_COMPONENTS] 
)
static

Compute the interpolated value of a vector variable.

Template Parameters
NUM_COMPONENTSNumber of components for the vector variable.

Compute the interpolated value of a variable.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
Parameters
NArray (for each support point) of shape function values at the coordinate the variable is to be interpolated.
varArray of variable values for each support point.
valueThe interpolated value of var.

This is the standard finite element interpolation operator of a discrete variable defined at the support points. The operator is expressed as:

\[ value = \sum_a^{numSupport} \left ( N_a var_a \right ), \]

Note
The shape function values N must be evaluated prior to calling this function.

◆ valueAndGradient()

template<int NUM_SUPPORT_POINTS, typename GRADIENT_TYPE >
static GEOS_HOST_DEVICE void geos::finiteElement::FiniteElementBase::valueAndGradient ( real64 const (&)  N[NUM_SUPPORT_POINTS],
GRADIENT_TYPE const &  gradN,
real64 const (&)  var[NUM_SUPPORT_POINTS],
real64 value,
real64(&)  gradVar[3] 
)
static

Calculate the value and gradient of a scalar valued support field at a point using the input basis function gradients.

Template Parameters
NUM_SUPPORT_POINTSThe number of support points for the element.
GRADIENT_TYPEThe type of the array object holding the shape
Parameters
NArray (for each support point) of shape function values at the coordinate the variable is to be interpolated.
gradNThe basis function gradients at a point in the element.
varThe vector valued support field that the gradient operator will be applied to.
valueThe value at the point for which N was specified.
gradVarThe gradient at the point for which gradN was specified.

Member Data Documentation

◆ m_viewDetJ

arrayView2d< real64 const > geos::finiteElement::FiniteElementBase::m_viewDetJ
protected

View to potentially hold pre-calculated weighted jacobian transformation determinants.

Definition at line 771 of file FiniteElementBase.hpp.


The documentation for this class was generated from the following file: