Function Operators

StandardFunctionOperators

StandardFunctionOperators are abstract types that represent fundamental (linear) operators, such as the identity, gradient, divergence, and others. These operators provide a unified interface for evaluating finite element basis functions in various ways.

List of Primitive Operators

StandardFunctionOperatorDescriptionMathematically
Identityidentity$v \rightarrow v$
IdentityComponent{c}identity of c-th component$v \rightarrow v_c$
NormalFluxnormal flux (function times normal)$v \rightarrow v \cdot \vec{n}$ (only ON_FACES)
TangentFluxtangent flux (function times tangent)$v \rightarrow v \cdot \vec{t}$ (only ON_EDGES)
Gradientgradient/Jacobian (as a vector)$v \rightarrow \nabla v$
SymmetricGradientsymmetric part of the gradient$v \rightarrow Voigt(\mathrm{sym}(\nabla v))$
Divergencedivergence$v \rightarrow \mathrm{div}(v) = \nabla \cdot v$
CurlScalarcurl operator 1D to 2D (rotated gradient)$v \rightarrow [-dv/dx_2,dv/dx_1]$
Curl2Dcurl operator 2D to 1D$v \rightarrow dv_1/dx_2 - dv_2/dx_1$
Curl3Dcurl operator 3D to 3D$v \rightarrow \nabla \times v$
HessianHesse matrix = all 2nd order derivatives (as a vector)$v \rightarrow D^2 v$ (e.g. in 2D: xx,xy,yx,yy for each component)
SymmetricHessian{a}symmetric part of Hesse matrix, offdiagonals scaled by a$v \rightarrow sym(D^2 v)$ (e.g. in 2D: xx,yy,a*xy for each component)
LaplacianLaplace Operator (diagonal of Hessian)$v \rightarrow \Delta v$ (e.g. in 2D: xx,yy for each component)
Note

The transformation from the reference domain to the physical domain differs for each finite element class. As a result, the evaluation of each function operator must be implemented specifically for every finite element class. Not all function operators are currently available for every dimension or element type, but new implementations are added as needed or upon request.

Additionally, function operators can be combined with user-defined kernels to postprocess/construct more advanced operators from the available primitives (for example, the deviatoric part of a tensor).

ExtendableFEMBase.Curl2DType
abstract type Curl2D <: ??

evaluates the curl of some two-dimensional vector field, i.e. Curl2D((u1,u2)) = du2/dx1 - du1/dx2

source
ExtendableFEMBase.Curl3DType
abstract type Curl3D <: ??

evaluates the curl of some three-dimensional vector field, i.e. Curl3D(u) = abla imes u

source
ExtendableFEMBase.NormalFluxType
abstract type NormalFlux <: ??

evaluates the normal-flux of the finite element function.

only available on FACES/BFACES and currently only for H1 and Hdiv elements

source
ExtendableFEMBase.OperatorPairType
abstract type OperatorPair{<:StandardFunctionOperator,<:StandardFunctionOperator} <: StandardFunctionOperator

allows to evaluate two operators in place of one, e.g. OperatorPair{Identity,Gradient}.

source
ExtendableFEMBase.OperatorTripleType
abstract type OperatorTriple{<:StandardFunctionOperator,<:StandardFunctionOperator} <: StandardFunctionOperator

allows to evaluate three operators in place of one, e.g. OperatorTriple{Identity,Gradient,Hessian}.

source
ExtendableFEMBase.SymmetricGradientType
abstract type SymmetricGradient{offdiagval} <: ??

evaluates the symmetric part of the gradient of the finite element function and returns its Voigt compression (off-diagonals on position [j,k] and [k,j] are added together and weighted with offdiagval).

source
ExtendableFEMBase.SymmetricHessianType
abstract type SymmetricHessian{offdiagval} <: ??

evaluates only the symmetric part of the Hessian of some (possibly vector-valued) finite element function. A concatenation of Voigt compressed Hessians for each component of the finite element functions is returned. The weighting of the mixed derivatives can be steered with offdiagval (e.g. √2 or 1 depending on the use case).

source
ExtendableFEMBase.Length4OperatorMethod
Length4Operator(
    _::Type{<:id},
    xdim::Int64,
    ncomponents::Int64
) -> Int64

number of expected components for the operator in xdim dimensions and a function value dimension of ncomponents

source

ReconstructionOperators

Special operators are provided to evaluate a primitive operator on a reconstructed version of a test function. These are useful for advanced discretizations and post-processing.

ExtendableFEMBase.ReconstructType
abstract type Reconstruct{FETypeR, O} <: ??

reconstruction operator: evaluates a reconstructed version of the finite element function.

FETypeR specifies the reconstruction space (needs to be defined for the finite element that it is applied to). O specifies the StandardFunctionOperator that shall be evaluated.

source

Divergence-Free Reconstruction Operators

For gradient-robust discretizations of certain classical non-divergence-conforming ansatz spaces, reconstruction operators are available that map a discretely divergence-free H1 function to a pointwise divergence-free H(div) function. Currently, such operators are implemented for the vector-valued Crouzeix-Raviart (H1CR) and Bernardi–Raugel (H1BR) finite element types, as well as for the P2-bubble (H1P2B) element in two dimensions.

Example: Reconst{HDIVRT0{d}, Identity} reconstructs the identity operator into HDIVRT0, and is available for H1BR{d} and H1CR{d} for d = 1, 2.

Operator Pairs (Experimental)

Two function operators can be combined into an OperatorPair, allowing you to provide two operators in each argument of an assembly pattern. Both operators must be well-defined on the relevant element geometries and finite element spaces, and their actions must be compatible with the input and result fields. This feature is experimental and may have limitations in some cases. An OperatorTriple is also available for combining three operators.