Function Operators

StandardFunctionOperators

StandardFunctionOperators are abstract types that encode primitive (linear) operators (like Identity, Gradient etc.) used to dispatch different evaluations of finite element basis functions.

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

As each finite element type is transformed differently from the reference domain to the general domain, the evaluation of each function operator has to be implemented for each finite element class. Currently, not every function operator works in any dimension and for any finite element. More evaluations are added as soon as they are needed (and possibly upon request). Also, the function operators can be combined with user-defined actions to evaluate other operators that can be build from the ones available (e.g. the deviator).

ExtendableFEMBase.Curl2DType
abstract type Curl2D <: StandardFunctionOperator

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

source
ExtendableFEMBase.Curl3DType
abstract type Curl3D <: StandardFunctionOperator

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

source
ExtendableFEMBase.CurlScalarType
abstract type CurlScalar <: StandardFunctionOperator

evaluates the curl of some scalar function in 2D, i.e. the rotated gradient.

source
ExtendableFEMBase.HessianType
abstract type Hessian <: StandardFunctionOperator

evaluates the full Hessian of some (possibly vector-valued) finite element function.

source
ExtendableFEMBase.LaplacianType
abstract type Laplacian <: StandardFunctionOperator

evaluates the Laplacian of some (possibly vector-valued) finite element function.

source
ExtendableFEMBase.NormalFluxType
abstract type NormalFlux <: StandardFunctionOperator

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} <: StandardFunctionOperator

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} <: StandardFunctionOperator

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

ReconstructionOperators

There are special operators that allow to evaluate a primitive operator of some discrete reconstructed version of a testfunction.

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

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 discretisations 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 Hdiv function. So far such operators are available for the vector-valued Crouzeix-Raviart (H1CR) and Bernardi–Raugel (H1BR) finite element types, as well as for the P2-bubble (H1P2B) finite element type in two dimensions.

Example: Reconst{HDIVRT0{d}, Identity} gives the reconstruction of 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 put into an OperatorPair so that one can provide effectively two operators in each argument of an assembly pattern. However, the user should make sure that both operators can be evaluated together reasonably (meaning both should be well-defined on the element geometries and the finite element space where the argument will be evaluated, and the action of the operator has to operate with coressponding input and result fields). This feature is still experimental and might have issues in some cases. OperatorTriple for a combination of three operators is also available.