GEOS
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
geos::finiteElement::BB_Tetrahedron< ORDER > Class Template Referencefinal

#include <BB_Tetrahedron.hpp>

Inheritance diagram for geos::finiteElement::BB_Tetrahedron< ORDER >:
Inheritance graph
[legend]

Public Member Functions

virtual ~BB_Tetrahedron ()=default
 
virtual GEOS_HOST_DEVICE localIndex getNumQuadraturePoints () const override
 Virtual getter for the number of quadrature points per element. More...
 
GEOS_HOST_DEVICE virtual GEOS_FORCE_INLINE localIndex getNumSupportPoints () const override
 Virtual getter for the number of support points per element. More...
 
GEOS_HOST_DEVICE virtual GEOS_FORCE_INLINE localIndex getMaxSupportPoints () const override
 Get the maximum number of support points for this element. More...
 
- Public Member Functions inherited from geos::finiteElement::FiniteElementBase
 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...
 
template<typename LEAF >
GEOS_HOST_DEVICE localIndex numSupportPoints (typename LEAF::StackVariables const &stack) const
 Getter for the number of support points per 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...
 
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...
 

Static Public Member Functions

GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumQuadraturePoints (StackVariables const &stack)
 Get the number of quadrature points. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex getNumSupportPoints (StackVariables const &stack)
 Get the number of support points. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 jacobianDeterminant (real64 const (&X)[4][3])
 Returns the determinant of the Jacobian of the element. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 faceJacobianDeterminant (localIndex face, real64 const (&X)[4][3])
 Calculate the determinant of the jacobian on the face opposite to the given vertex. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN (localIndex const, real64(&N)[numNodes])
 Calculate shape functions values for each support point at a quadrature point. More...
 
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 point. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN (real64 const (&lambda)[4], real64(&N)[numNodes])
 Calculate shape functions values at a single point using De Casteljau's algorithm. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcN (real64 const (&X)[4][3], real64 const (&coords)[3], real64(&N)[numNodes])
 Calculate shape functions values at a single point, given the coordinates of the tetrahedron vertices, using De Casteljau's algorithm. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcNandGradN (real64 const (&lambda)[4], real64 const (&N)[numNodes], real64(&gradN)[numNodes][4])
 Calculate the values and derivatives of shape functions with respect to barycentric coordinates at a single point using De Casteljau's algorithm. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void calcNandGradN (real64 const (&X)[4][3], real64 const (&coords)[3], real64(&N)[numNodes], real64(&gradN)[numNodes][4])
 Calculate the shape functions values and derivatives at a single point, given the coorginates of the tetrahedron vertices, using De Casteljau's algorithm. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 calcGradN (localIndex const q, real64 const (&X)[numNodes][3], real64(&gradN)[numNodes][3])
 Calculate the shape functions derivatives wrt the physical coordinates. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 calcGradN (localIndex const q, real64 const (&X)[numNodes][3], StackVariables const &stack, real64(&gradN)[numNodes][3])
 Calculate the shape functions derivatives wrt the physical coordinates. More...
 
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 integralTerm (const int a, const int b, const int c)
 Computes a! / ( b! * c! ) with b + c >= a >= b, c. More...
 
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 correctionFactorDerivative (int const i1, int const j1, int const k1, int const l1, int const i2, int const j2, int const k2, int const l2, int const dim)
 Computes the correction factor for the superposition integral in case a function has been derived once. The indices of the original function are ii1, j1, k1 and l1, and those of the once-derived one are i2, j2, k2 and l2. More...
 
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 computeSuperpositionIntegral (const int i1, const int j1, const int k1, const int l1, const int i2, const int j2, const int k2, const int l2)
 Computes the superposition integral between Bernstein-Bézier functions indexed by (i1, j1, k1, l1) and (i1, j2, k2, l2) More...
 
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 computeFaceSuperpositionIntegral (const int i1, const int j1, const int k1, const int i2, const int j2, const int k2)
 Computes the superposition integral over a face between Bernstein-Bézier functions whose indices are given by (i1, j1, k1, l1=0 ) and (i1, j2, k2, l2=0) More...
 
template<int I, int J, int K>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int dofIndex ()
 Computes the local degree of freedom index given the shape function indices (i, j, k, l) for each vertex. Only i, j and k are needed since i + j + k + l = order. More...
 
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int dofIndex (const int i, const int j, const int k)
 Computes the local degree of freedom index given the shape function indices (i, j, k, l) for each vertex. Only i, j and k are needed since i + j + k + l = order This version takes indices as parameters and can be called with non-constexpr index values. More...
 
template<int C, int VTX>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int indexToIJKL ()
 Computes the local degree of freedom index given the shape function indices for each vertex. More...
 
template<typename FUNC , int... Is>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void loop_impl (FUNC &&func, std::integer_sequence< int, Is... >)
 Helper function for static for loop. More...
 
template<int N, typename FUNC >
static constexpr GEOS_HOST_DEVICE void loop (FUNC const &func)
 Helper function for static for loop. More...
 
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void barycentricCoordinateLoop (FUNC &&func)
 Helper function for loop over barycentric coordinates. More...
 
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void basisLoop (FUNC &&func)
 Helper function for loop over tet basis functions. More...
 
template<int c1, int i1, int j1, int k1, int l1, typename F , int... Is>
static constexpr GEOS_HOST_DEVICE void call_matching_cases (F &&func, std::integer_sequence< int, Is... >)
 Helper function for loop over tet basis functions that have one index in a given set of indices. If multiple indices are in the given list, the callback is called multiple times. This herlper is useful to avoid too much fold operations. More...
 
template<int... Is, typename FUNC >
static constexpr GEOS_HOST_DEVICE void conditionalBasisLoop (FUNC const &func)
 Helper function for loop over tet basis functions that have one index in a given set of indices. If multiple indices are in the given list, the callback is called multiple times. More...
 
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void faceBarycentricCoordinateLoop (FUNC &&func)
 Helper function for loop over barycentric coordinates of a face. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void computeReferenceMassMatrix (arraySlice2d< real64 > const &m)
 Computes the reference mass matrix, i.e., the superposition matrix of the shape functions in barycentric coordinates. The real-world mass matrix can be obtained by using the multiplying this matrix by the determinant of the Jacobian. More...
 
template<typename DAMPING >
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void computeReferenceDampingMatrix (real64(&d)[numNodes][numNodes], bool const face1Damped, bool const face2Damped, bool const face3Damped, bool const face4Damped)
 Computes the reference damping matrix, i.e., the superposition matrix of the shape functions in barycentric coordinates over faces. The real-world mass matrix can be obtained by using the multiplying this matrix by the determinant of the face Jacobian. More...
 
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void computeMassTerm (real64 const (&X)[4][3], FUNC &&func)
 Computes the non-zero contributions inside the element of the mass matrix M, i.e., the superposition matrix of shape functions. More...
 
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 computeFluxDerivativeFactor (real64 const (&X)[4][3], int x1, int x2, int o1, int o2)
 Function to compute the factor for the flux derivative term. More...
 
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void computeStiffnessTerm (real64 const (&X)[4][3], FUNC &&func)
 Computes the non-zero contributions inside the element of the stiffness matrix R, i.e., the superposition matrix of first derivatives of shape functions. More...
 
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 edgeLength2 (localIndex i1, localIndex i2, real64 const (&X)[4][3])
 
template<typename FUNCP , typename FUNCF >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void computeSurfaceTerms (real64 const (&X)[4][3], FUNCP &&funcP, FUNCF &&funcF)
 Computes the non-zero contributions inside the element of the surface terms, including the value of the superposition integral of basis functions (used for the penalty and damping terms) and the superposition integral of the derivative of a function with the value of the other, used for the flux terms. More...
 
- Static Public Member Functions inherited from geos::finiteElement::FiniteElementBase
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...
 
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...
 
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...
 
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...
 
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...
 

Static Public Attributes

static constexpr int order = ORDER
 The order of the finite element.
 
constexpr static localIndex numNodes = ( ORDER + 1 ) * ( ORDER + 2 ) * ( ORDER + 3 ) / 6
 The number of shape functions per element.
 
constexpr static localIndex numNodesPerFace = ( ORDER + 1 ) * ( ORDER + 2 ) / 2
 The number of shape functions per face.
 
constexpr static localIndex maxSupportPoints = numNodes
 The maximum number of support points per element.
 
constexpr static localIndex numQuadraturePoints = numNodes
 The number of quadrature points per element.
 
- Static Public Attributes inherited from geos::finiteElement::FiniteElementBase
constexpr static int numSamplingPointsPerDirection = 10
 Number of sampling points.
 

Additional Inherited Members

- Protected Attributes inherited from geos::finiteElement::FiniteElementBase
arrayView4d< real64 const > m_viewGradN
 View to potentially hold pre-calculated shape function gradients.
 
arrayView2d< real64 const > m_viewDetJ
 

Detailed Description

template<int ORDER>
class geos::finiteElement::BB_Tetrahedron< ORDER >

This class contains the kernel accessible functions specific to the Bernstein-Bézier (BB) modal any-order tetrahedron finite element with Gaussian quadrature rules. Available functions are tailored for Discontinuous Galerkin (DG) applications. In barycentric coordinates l1, l2, l3, l4, the function indexed by (i1, i2, i3, i4) corresponds to the function (i1+i2+i3+i4+3)! / (i1! i2! i3! i4!) l1^i1 l2^i2 l3^i3 l4^i4 and they integrate to one over the reference element defined by 0<=l1,l2,l3,l3<=1, l1+l2+l3+l4=1

Definition at line 45 of file BB_Tetrahedron.hpp.

Constructor & Destructor Documentation

◆ ~BB_Tetrahedron()

template<int ORDER>
virtual geos::finiteElement::BB_Tetrahedron< ORDER >::~BB_Tetrahedron ( )
virtualdefault

Doxygen_Suppress

Member Function Documentation

◆ barycentricCoordinateLoop()

template<int ORDER>
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::barycentricCoordinateLoop ( FUNC &&  func)
inlinestaticconstexpr

Helper function for loop over barycentric coordinates.

Template Parameters
FUNCthe callback function
Parameters
functhe callback function to call for each index

Definition at line 678 of file BB_Tetrahedron.hpp.

◆ basisLoop()

template<int ORDER>
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::basisLoop ( FUNC &&  func)
inlinestaticconstexpr

Helper function for loop over tet basis functions.

Template Parameters
FUNCthe callback function
Parameters
functhe callback function to call for each index

Definition at line 691 of file BB_Tetrahedron.hpp.

◆ calcGradN() [1/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::calcGradN ( localIndex const  q,
real64 const (&)  X[numNodes][3],
real64(&)  gradN[numNodes][3] 
)
inlinestatic

Calculate the shape functions derivatives wrt the physical coordinates.

Parameters
qIndex of the quadrature point.
XArray containing the coordinates of the mesh support points.
gradNArray to contain the shape function derivatives for all support points at the coordinates of the quadrature point q.
Returns
The determinant of the parent/physical transformation matrix.

Definition at line 411 of file BB_Tetrahedron.hpp.

◆ calcGradN() [2/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::calcGradN ( localIndex const  q,
real64 const (&)  X[numNodes][3],
StackVariables const &  stack,
real64(&)  gradN[numNodes][3] 
)
inlinestatic

Calculate the shape functions derivatives wrt the physical coordinates.

Parameters
qIndex of the quadrature point.
XArray containing the coordinates of the mesh support points.
stackVariables allocated on the stack as filled by setupStack.
gradNArray to contain the shape function derivatives for all support points at the coordinates of the quadrature point q.
Returns
The determinant of the parent/physical transformation matrix.

Definition at line 434 of file BB_Tetrahedron.hpp.

◆ calcN() [1/4]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcN ( localIndex const  q,
StackVariables const &  stack,
real64(&)  N[numNodes] 
)
inlinestatic

Calculate shape functions values for each support point at a quadrature point.

Parameters
qIndex of the quadrature point.
stackVariables allocated on the stack as filled by setupStack.
NAn array to pass back the shape function values for each support point.

Definition at line 194 of file BB_Tetrahedron.hpp.

◆ calcN() [2/4]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcN ( localIndex const  ,
real64(&)  N[numNodes] 
)
inlinestatic

Calculate shape functions values for each support point at a quadrature point.

Parameters
NAn array to pass back the shape function values for each support point.

Definition at line 174 of file BB_Tetrahedron.hpp.

◆ calcN() [3/4]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcN ( real64 const (&)  lambda[4],
real64(&)  N[numNodes] 
)
inlinestatic

Calculate shape functions values at a single point using De Casteljau's algorithm.

Parameters
[in]lambdabarycentric coordinates of the point in thetetrahedron
[out]NThe shape function.

Definition at line 208 of file BB_Tetrahedron.hpp.

◆ calcN() [4/4]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcN ( real64 const (&)  X[4][3],
real64 const (&)  coords[3],
real64(&)  N[numNodes] 
)
inlinestatic

Calculate shape functions values at a single point, given the coordinates of the tetrahedron vertices, using De Casteljau's algorithm.

Parameters
[in]XArray containing the coordinates of the support points.
[in]coordsThe parent coordinates at which to evaluate the shape function value, in the reference element
[in]NAn array to pass back the shape function values for each support point

Definition at line 256 of file BB_Tetrahedron.hpp.

◆ calcNandGradN() [1/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcNandGradN ( real64 const (&)  lambda[4],
real64 const (&)  N[numNodes],
real64(&)  gradN[numNodes][4] 
)
inlinestatic

Calculate the values and derivatives of shape functions with respect to barycentric coordinates at a single point using De Casteljau's algorithm.

Parameters
[in]lambdabarycentric coordinates of the point in thetetrahedron
[out]NThe shape function values.
[out]gradNThe derivatives of the shape functions with respect to the lambdas

Definition at line 295 of file BB_Tetrahedron.hpp.

◆ calcNandGradN() [2/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::calcNandGradN ( real64 const (&)  X[4][3],
real64 const (&)  coords[3],
real64(&)  N[numNodes],
real64(&)  gradN[numNodes][4] 
)
inlinestatic

Calculate the shape functions values and derivatives at a single point, given the coorginates of the tetrahedron vertices, using De Casteljau's algorithm.

Parameters
[in]XAn array containing the coordinates of the tetrahedra
[in]coordsThe parent coordinates at which to evaluate the shape function value, in the reference element
[in]NAn array to pass back the shape function values for each support
[out]gradNThe derivatives of the shape functions with respect to the lambdas

Definition at line 354 of file BB_Tetrahedron.hpp.

◆ call_matching_cases()

template<int ORDER>
template<int c1, int i1, int j1, int k1, int l1, typename F , int... Is>
static constexpr GEOS_HOST_DEVICE void geos::finiteElement::BB_Tetrahedron< ORDER >::call_matching_cases ( F &&  func,
std::integer_sequence< int, Is... >   
)
inlinestaticconstexpr

Helper function for loop over tet basis functions that have one index in a given set of indices. If multiple indices are in the given list, the callback is called multiple times. This herlper is useful to avoid too much fold operations.

Template Parameters
c1the dof index in the element
i1the index with respect to the first vertex
j1the index with respect to the second vertex
k1the index with respect to the third vertex
l1the index with respect to the fourth vertex
Fthe callback function type
Isthe set of indices to check against
Parameters
functhe callback function to call for each matching index

Definition at line 739 of file BB_Tetrahedron.hpp.

◆ computeFaceSuperpositionIntegral()

template<int ORDER>
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::computeFaceSuperpositionIntegral ( const int  i1,
const int  j1,
const int  k1,
const int  i2,
const int  j2,
const int  k2 
)
inlinestaticconstexpr

Computes the superposition integral over a face between Bernstein-Bézier functions whose indices are given by (i1, j1, k1, l1=0 ) and (i1, j2, k2, l2=0)

Parameters
[in]i1
[in]j1
[in]k1
[in]i2
[in]j2
[in]k2
Returns
the superposition integral over the barycentric coordinates

Definition at line 536 of file BB_Tetrahedron.hpp.

◆ computeFluxDerivativeFactor()

template<int ORDER>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::computeFluxDerivativeFactor ( real64 const (&)  X[4][3],
int  x1,
int  x2,
int  o1,
int  o2 
)
inlinestaticconstexpr

Function to compute the factor for the flux derivative term.

Parameters
XArray containing the coordinates of the support points.
x1Index of the first edge vertex
x2Index of the second edge vertex
o1Index of the first face vertex
o2Index of the second face vertex
Returns
The factor for the flux derivative term

Definition at line 1127 of file BB_Tetrahedron.hpp.

◆ computeMassTerm()

template<int ORDER>
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::computeMassTerm ( real64 const (&)  X[4][3],
FUNC &&  func 
)
inlinestaticconstexpr

Computes the non-zero contributions inside the element of the mass matrix M, i.e., the superposition matrix of shape functions.

Parameters
XArray containing the coordinates of the support points.
funcCallback function accepting three parameters: i, j (local d.o.f. inside the element) and M_ij

Definition at line 1068 of file BB_Tetrahedron.hpp.

◆ computeReferenceDampingMatrix()

template<int ORDER>
template<typename DAMPING >
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::computeReferenceDampingMatrix ( real64(&)  d[numNodes][numNodes],
bool const  face1Damped,
bool const  face2Damped,
bool const  face3Damped,
bool const  face4Damped 
)
inlinestatic

Computes the reference damping matrix, i.e., the superposition matrix of the shape functions in barycentric coordinates over faces. The real-world mass matrix can be obtained by using the multiplying this matrix by the determinant of the face Jacobian.

Parameters
[out]dThe damping matrix
[in]face1DampedWhether the first face contributes to the damping term
[in]face2DampedWhether the second face contributes to the damping term
[in]face3DampedWhether the third face contributes to the damping term
[in]face4DampedWhether the fourth face contributes to the damping term

Definition at line 1027 of file BB_Tetrahedron.hpp.

◆ computeReferenceMassMatrix()

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::computeReferenceMassMatrix ( arraySlice2d< real64 > const &  m)
inlinestatic

Computes the reference mass matrix, i.e., the superposition matrix of the shape functions in barycentric coordinates. The real-world mass matrix can be obtained by using the multiplying this matrix by the determinant of the Jacobian.

Parameters
[out]mThe mass matrix

Definition at line 987 of file BB_Tetrahedron.hpp.

◆ computeStiffnessTerm()

template<int ORDER>
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::computeStiffnessTerm ( real64 const (&)  X[4][3],
FUNC &&  func 
)
inlinestaticconstexpr

Computes the non-zero contributions inside the element of the stiffness matrix R, i.e., the superposition matrix of first derivatives of shape functions.

Parameters
XArray containing the coordinates of the support points.
funcCallback function accepting three parameters: i, j (local d.o.f. inside the element) and R_ij

Definition at line 1175 of file BB_Tetrahedron.hpp.

◆ computeSuperpositionIntegral()

template<int ORDER>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::computeSuperpositionIntegral ( const int  i1,
const int  j1,
const int  k1,
const int  l1,
const int  i2,
const int  j2,
const int  k2,
const int  l2 
)
inlinestaticconstexpr

Computes the superposition integral between Bernstein-Bézier functions indexed by (i1, j1, k1, l1) and (i1, j2, k2, l2)

Parameters
[in]i1
[in]j1
[in]k1
[in]l1
[in]i2
[in]j2
[in]k2
[in]l2
Returns
the superposition integral over the barycentric coordinates

Definition at line 512 of file BB_Tetrahedron.hpp.

◆ computeSurfaceTerms()

template<int ORDER>
template<typename FUNCP , typename FUNCF >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::computeSurfaceTerms ( real64 const (&)  X[4][3],
FUNCP &&  funcP,
FUNCF &&  funcF 
)
inlinestaticconstexpr

Computes the non-zero contributions inside the element of the surface terms, including the value of the superposition integral of basis functions (used for the penalty and damping terms) and the superposition integral of the derivative of a function with the value of the other, used for the flux terms.

Parameters
XArray containing the coordinates of the support points.
funcPCallback function for non-zero penalty-type terms, accepting seven parameters: c1, c2 (local d.o.f. inside the element), f (index of the face, i.e., index of the opposite vertex for this element), i1, j1 and k1 (local indices for the first shape function) and value i2, j2 and k2 (local indices for the second shape function) and value
funcFCallback function for non-zero flux-type terms, accepting four parameters: c2, c2 (local d.o.f. inside the element), f (index of the face, i.e., index of the opposite vertex for this element), and value i1, j1 and k1 (local indices for the first shape function) and value i2, j2 and k2 (local indices for the second shape function) and value

Definition at line 1277 of file BB_Tetrahedron.hpp.

◆ conditionalBasisLoop()

template<int ORDER>
template<int... Is, typename FUNC >
static constexpr GEOS_HOST_DEVICE void geos::finiteElement::BB_Tetrahedron< ORDER >::conditionalBasisLoop ( FUNC const &  func)
inlinestaticconstexpr

Helper function for loop over tet basis functions that have one index in a given set of indices. If multiple indices are in the given list, the callback is called multiple times.

Template Parameters
FUNCthe callback function
Isthe setindices
Parameters
functhe callback function to call for each index

Definition at line 801 of file BB_Tetrahedron.hpp.

◆ correctionFactorDerivative()

template<int ORDER>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::correctionFactorDerivative ( int const  i1,
int const  j1,
int const  k1,
int const  l1,
int const  i2,
int const  j2,
int const  k2,
int const  l2,
int const  dim 
)
inlinestaticconstexpr

Computes the correction factor for the superposition integral in case a function has been derived once. The indices of the original function are ii1, j1, k1 and l1, and those of the once-derived one are i2, j2, k2 and l2.

Parameters
[in]i1
[in]j1
[in]k1
[in]l1
[in]i2
[in]j2
[in]k2
[in]l2
[in]dimthe dimension, 2 or 3
Returns
the correction factor to be applied to the superposition integral

Definition at line 491 of file BB_Tetrahedron.hpp.

◆ dofIndex() [1/2]

template<int ORDER>
template<int I, int J, int K>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int geos::finiteElement::BB_Tetrahedron< ORDER >::dofIndex ( )
inlinestaticconstexpr

Computes the local degree of freedom index given the shape function indices (i, j, k, l) for each vertex. Only i, j and k are needed since i + j + k + l = order.

Template Parameters
IThe index with respect to the first vertex
JThe index with respect to the second vertex
KThe index with respect to the third vertex
Returns
The local degree of freedom index

Definition at line 561 of file BB_Tetrahedron.hpp.

◆ dofIndex() [2/2]

template<int ORDER>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int geos::finiteElement::BB_Tetrahedron< ORDER >::dofIndex ( const int  i,
const int  j,
const int  k 
)
inlinestaticconstexpr

Computes the local degree of freedom index given the shape function indices (i, j, k, l) for each vertex. Only i, j and k are needed since i + j + k + l = order This version takes indices as parameters and can be called with non-constexpr index values.

Parameters
iThe index with respect to the first vertex
jThe index with respect to the second vertex
kThe index with respect to the third vertex
Returns
The local degree of freedom index

Definition at line 582 of file BB_Tetrahedron.hpp.

◆ edgeLength2()

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::edgeLength2 ( localIndex  i1,
localIndex  i2,
real64 const (&)  X[4][3] 
)
inlinestatic

Compute th length of the edge between two vertices i1 and i2

Parameters
i1Index of the first vertex
i2Index of the second vertex
XArray containing the coordinates of the support points.
Returns
The squared length of the edge between the two vertices

Definition at line 1247 of file BB_Tetrahedron.hpp.

◆ faceBarycentricCoordinateLoop()

template<int ORDER>
template<typename FUNC >
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::faceBarycentricCoordinateLoop ( FUNC &&  func)
inlinestaticconstexpr

Helper function for loop over barycentric coordinates of a face.

Template Parameters
FUNCthe callback function
Parameters
functhe callback function to call for each index

Definition at line 971 of file BB_Tetrahedron.hpp.

◆ faceJacobianDeterminant()

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::faceJacobianDeterminant ( localIndex  face,
real64 const (&)  X[4][3] 
)
inlinestatic

Calculate the determinant of the jacobian on the face opposite to the given vertex.

Parameters
[in]faceThe index of the vertex opposite to the desired face
[in]XThe coordinates of the tetrahedron
Returns
the (absolute value of the) determinant of the Jacobian on the face

Definition at line 146 of file BB_Tetrahedron.hpp.

◆ getMaxSupportPoints()

template<int ORDER>
GEOS_HOST_DEVICE virtual GEOS_FORCE_INLINE localIndex geos::finiteElement::BB_Tetrahedron< ORDER >::getMaxSupportPoints ( ) const
inlineoverridevirtual

Get the maximum number of support points for this element.

This should be used to know the size of pre-allocated objects whose size depend on the number of support points.

Returns
The number of maximum support points for this element.

Implements geos::finiteElement::FiniteElementBase.

Definition at line 98 of file BB_Tetrahedron.hpp.

◆ getNumQuadraturePoints() [1/2]

template<int ORDER>
virtual GEOS_HOST_DEVICE localIndex geos::finiteElement::BB_Tetrahedron< ORDER >::getNumQuadraturePoints ( ) const
inlineoverridevirtual

Virtual getter for the number of quadrature points per element.

Returns
The number of quadrature points per element.

Implements geos::finiteElement::FiniteElementBase.

Definition at line 71 of file BB_Tetrahedron.hpp.

◆ getNumQuadraturePoints() [2/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex geos::finiteElement::BB_Tetrahedron< ORDER >::getNumQuadraturePoints ( StackVariables const &  stack)
inlinestatic

Get the number of quadrature points.

Parameters
stackStack variables as filled by setupStack.
Returns
The number of quadrature points.

Definition at line 83 of file BB_Tetrahedron.hpp.

◆ getNumSupportPoints() [1/2]

template<int ORDER>
GEOS_HOST_DEVICE virtual GEOS_FORCE_INLINE localIndex geos::finiteElement::BB_Tetrahedron< ORDER >::getNumSupportPoints ( ) const
inlineoverridevirtual

Virtual getter for the number of support points per element.

Returns
The number of support points per element.

Implements geos::finiteElement::FiniteElementBase.

Definition at line 91 of file BB_Tetrahedron.hpp.

◆ getNumSupportPoints() [2/2]

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE localIndex geos::finiteElement::BB_Tetrahedron< ORDER >::getNumSupportPoints ( StackVariables const &  stack)
inlinestatic

Get the number of support points.

Parameters
stackObject that holds stack variables.
Returns
The number of support points.

Definition at line 110 of file BB_Tetrahedron.hpp.

◆ indexToIJKL()

template<int ORDER>
template<int C, int VTX>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE int geos::finiteElement::BB_Tetrahedron< ORDER >::indexToIJKL ( )
inlinestaticconstexpr

Computes the local degree of freedom index given the shape function indices for each vertex.

Template Parameters
CThe dof index in the element
VTXthe vertex with respect to
Returns
i, j, k, l if VTX= 0, 1, 2 or 3 resepctively (with i + j + k + l = order)

Definition at line 603 of file BB_Tetrahedron.hpp.

◆ integralTerm()

template<int ORDER>
GEOS_HOST_DEVICE constexpr static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::integralTerm ( const int  a,
const int  b,
const int  c 
)
inlinestaticconstexpr

Computes a! / ( b! * c! ) with b + c >= a >= b, c.

Parameters
[in]a
[in]b
[in]c
Returns
a!/(b!*c!)

Definition at line 452 of file BB_Tetrahedron.hpp.

◆ jacobianDeterminant()

template<int ORDER>
GEOS_HOST_DEVICE static GEOS_FORCE_INLINE real64 geos::finiteElement::BB_Tetrahedron< ORDER >::jacobianDeterminant ( real64 const (&)  X[4][3])
inlinestatic

Returns the determinant of the Jacobian of the element.

Parameters
[in]XThe coordinates of the tetrahedron
Returns
the (absolute value of the) determinant of the Jacobian on the element

Definition at line 124 of file BB_Tetrahedron.hpp.

◆ loop()

template<int ORDER>
template<int N, typename FUNC >
static constexpr GEOS_HOST_DEVICE void geos::finiteElement::BB_Tetrahedron< ORDER >::loop ( FUNC const &  func)
inlinestaticconstexpr

Helper function for static for loop.

Template Parameters
Nthe number of iterations
FUNCthe callback function
Parameters
functhe callback function to call for each index This function recursively calls itself until N reaches 0, at which point it stops.

Definition at line 659 of file BB_Tetrahedron.hpp.

◆ loop_impl()

template<int ORDER>
template<typename FUNC , int... Is>
GEOS_HOST_DEVICE static constexpr GEOS_FORCE_INLINE void geos::finiteElement::BB_Tetrahedron< ORDER >::loop_impl ( FUNC &&  func,
std::integer_sequence< int, Is... >   
)
inlinestaticconstexpr

Helper function for static for loop.

Template Parameters
FUNCthe callback function
...Isinteger indices of the loop
Parameters
functhe callback function to call for each index
Template Parameters
Isthe integer indices of the loop

Definition at line 645 of file BB_Tetrahedron.hpp.


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