FESpace

A finite element space (FESpace) represents the set of all functions that can be expressed using a particular finite element type on a given computational grid. In ExtendableFEMBase, constructing a finite element space requires only specifying the finite element type and the grid; all necessary degree-of-freedom (dof) mappings are generated automatically on first access.

ExtendableFEMBase.FESpaceType
struct FESpace{Tv, Ti, FEType<:AbstractFiniteElement,AT<:AssemblyType}
	name::String                          
	broken::Bool                         
	ndofs::Int                            
	coffset::Int                          
	xgrid::ExtendableGrid[Tv,Ti}           
	dofgrid::ExtendableGrid{Tv,Ti}	      
	dofmaps::Dict{Type{<:AbstractGridComponent},Any} 
    interpolators::Dict{Type{<:AssemblyType}, Any} 
end

A finite element space representing the global collection of degrees of freedom (DOFs) for a given finite element type and assembly type on a computational grid.

FESpace encapsulates the mapping between mesh entities (cells, faces, edges, etc.) and DOFs, as well as metadata and auxiliary structures needed for assembly, interpolation, and evaluation of finite element functions.

Type Parameters

  • Tv: Value type for grid coordinates (e.g., Float64).
  • Ti: Integer type for grid indices (e.g., Int64).
  • FEType: The finite element type (e.g., H1P1).
  • AT: The assembly type (e.g., ON_CELLS).

Fields

  • name::String: Name of the finite element space.
  • broken::Bool: Whether the space is "broken" (discontinuous across elements).
  • ndofs::Int64: Total number of global degrees of freedom.
  • coffset::Int: Offset for component DOFs (for vector-valued or mixed spaces).
  • xgrid::ExtendableGrid{Tv, Ti}: Reference to the computational grid.
  • dofgrid::ExtendableGrid{Tv, Ti}: Grid used for DOF numbering (may be a subgrid of xgrid).
  • dofmaps::Dict{Type{<:AbstractGridComponent}, Any}: Dictionary of DOF maps for different grid components (cells, faces, edges, etc.).
  • interpolators::Dict{Type{<:AssemblyType}, Any}: Dictionary of interpolation operators for different assembly types.
source
ExtendableFEMBase.FESpaceMethod
FESpace{FEType, AT}(xgrid::ExtendableGrid{Tv, Ti}; name = "", regions = nothing, broken::Bool = false)

Constructs a finite element space (FESpace) of the given finite element type (FEType) and assembly type (AT) on the provided computational grid.

  • The FESpace represents the global collection of degrees of freedom (DOFs) for the specified finite element type and assembly type, and manages the mapping between mesh entities (cells, faces, edges, etc.) and DOFs.
  • The broken switch allows the creation of a "broken" (discontinuous) finite element space, where basis functions are not required to be continuous across element boundaries.
  • The regions argument can be used to restrict the space to a subset of the grid.
  • If no AT is provided, the space is generated on cells (ON_CELLS).

Arguments

  • xgrid::ExtendableGrid{Tv, Ti}: The computational grid on which the finite element space is defined.

Keyword Arguments

  • `name::String: Name for the finite element space (for identification and debugging) (default: "", triggers automatic generation from FEType and broken).
  • regions: Optional subset of the grid to restrict the space (default: all regions).
  • broken: Whether to create a broken (discontinuous) space (default: false).
source
Base.eltypeMethod
eltype(
    _::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT}
) -> Type{FEType} where FEType<:AbstractFiniteElement

Custom eltype function for FESpace returns the finite element type parameter of the finite element space.

source
Base.get!Method
get!(FES::FESpace, DM::Type{<:DofMap}) -> Any

To be called by getindex. This triggers lazy creation of non-existing dofmaps

source
Base.getindexMethod
Base.getindex(FES::FESpace,DM::Type{<:DofMap})

Generic method for obtaining dofmap. This method is mutating in the sense that non-existing dofmaps are created on demand. Due to the fact that components are stored as Any the return value triggers type instability.

source
Base.showMethod
show(
    io::IO,
    FES::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT}
)

Custom show function for FESpace that prints some information and all available dofmaps.

source
ExtendableFEMBase.assemblytypeMethod
assemblytype(
    _::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT}
) -> Any

returns the assembly type parameter of the finite element space, i.e. on which entities of the grid the finite element is defined.

source
ExtendableFEMBase.coffsetMethod
coffset(FES::FESpace) -> Int64

returns the offset between the degrees of freedom of each component (i.e. the number of scalar degrees of freedom that influence a component, vector-valued degrees of freedom are stored at the end).

source
ExtendableFEMBase.count_ndofsMethod
count_ndofs(
    xgrid,
    FEType,
    broken::Bool
) -> Tuple{Int64, Int64}

counts the total number of degrees of freedom for the FEType for the whole grid

source
ExtendableFEMBase.get_basisMethod
get_basis(
    AT::Type{<:ExtendableGrids.AssemblyType},
    FEType::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#233"

returns the closure function of form

closure(refbasis, xref)

that computes the evaluations of the basis functions on the reference geometry

source
ExtendableFEMBase.get_basissubsetMethod
get_basissubset(
    ::Type{<:ExtendableGrids.AssemblyType},
    FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    ::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#261"{_A, Int64} where _A
get_basissubset(
    ::Type{<:ExtendableGrids.AssemblyType},
    FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    ::Type{<:ExtendableGrids.AbstractElementGeometry},
    xgrid
) -> Union{ExtendableFEMBase.var"#21#22", ExtendableFEMBase.var"#closure#270"{_A, _B, _C, Int64, Int64} where {_A, _B, _C}}

returns a closure function of the form

closure(subset_ids::Array{Int, 1}, cell)

which returns the ids of the local basis functions needed on the cell. Usually, subset_ids = 1:ndofs (meaning that all basis functions on the reference cells are used),

See e.g. the 3D implementation of BDM1 or H1P3 where different basis functions are chosen depending on the face orientations (which in 3D is not just a sign)

source
ExtendableFEMBase.get_coefficientsMethod
get_coefficients(
    ::Type{<:ExtendableGrids.AssemblyType},
    FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    ::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#221"
get_coefficients(
    ::Type{<:ExtendableGrids.AssemblyType},
    FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    ::Type{<:ExtendableGrids.AbstractElementGeometry},
    xgrid
) -> ExtendableFEMBase.var"#closure#287"

returns the coefficients for local evaluations of finite element functions ( see e.g. h1v_br.jl for a use-case)

source
ExtendableFEMBase.get_ndofsMethod
get_ndofs(
    AT::Type{<:ExtendableGrids.AssemblyType},
    FEType::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Any

returns the number of degrees of freedom for the given AssemblyType, FEType and geometry

source
ExtendableFEMBase.get_ndofs_allMethod
get_ndofs_all(
    AT::Type{<:ExtendableGrids.AssemblyType},
    FEType::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Int64

returns the total number of degrees of freedom for the given AssemblyType, FEType and geometry

source
ExtendableFEMBase.get_polynomialorderMethod
get_polynomialorder(
    AT::Type{<:ExtendableGrids.AssemblyType},
    FE::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Int64

returns the expected polynomial order of the basis functions (used to determine default quadrature orders)

source
ExtendableFEMBase.interior_dofs_offsetMethod
interior_dofs_offset(
    AT::Type{<:ExtendableGrids.AssemblyType},
    FE::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Int64

returns the offset for the interior degrees of freedom (all vertex, face and edge dofs are assumed to come first and before that number)

source
ExtendableFEMBase.isdefinedMethod
isdefined(
    FEType::Type{<:AbstractFiniteElement},
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Bool

tells if FEType is defined on this ElementGeometry

source
ExtendableFEMBase.xgridMethod
xgrid(FES::FESpace) -> ExtendableGrids.ExtendableGrid

returns the computational grid of the finite element space.

source

DofMaps

ExtendableFEMBase.DofMapType
abstract type DofMap <: ExtendableGrids.AbstractGridAdjacency

A DofMap encodes the association of global DoF indices to mesh entities of a given type (e.g., cells, faces, edges) for a particular finite element space. This mapping is essential for assembling system matrices and vectors, applying boundary conditions, and extracting solution values. DofMaps are typically not constructed directly by users. Instead, they are generated/managed internally by the finite element space and accessed as needed for assembly or evaluation tasks.

Implementation

  • All concrete DofMap subtypes (e.g., CellDofs, FaceDofs, EdgeDofs, etc.) specify the mesh entity type to which DoFs are attached.
  • DofMaps are stored as ExtendableGrids.AbstractGridAdjacency (usually VariableTargetAdjacency, SerialVariableTargetAdjacency or AbstractGridIntegerArray2D) objects within the finite element space and are generated automatically as needed.
  • The appropriate DofMap for a given assembly type can be accessed via FESpace[DofMapSubtype].

Example

cell_dofs = FESpace[CellDofs]   # Get the cell-based DoF map
face_dofs = FESpace[FaceDofs]   # Get the face-based DoF map
source
ExtendableFEMBase.ParsedDofMapMethod
ParsedDofMap(
    pattern::String,
    ncomponents,
    EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.ParsedDofMap

parses a dofmap string (defined in each finite element definition file) and generated a ParsedDofMap for a certain number of components and the given geometry

source
ExtendableFEMBase.Dofmap4AssemblyTypeMethod
Dofmap4AssemblyType(
    FES::FESpace,
    AT::Type{<:ExtendableGrids.AssemblyType}
) -> Any

yields the coressponding dofmap of the FESpace for a given AssemblyType (assumed with respect to the (parent)grid of the FESpace)

source
ExtendableFEMBase.EffAT4AssemblyTypeMethod
EffAT4AssemblyType(
    _::Type{ON_CELLS},
    _::Type{ON_CELLS}
) -> Type{ON_CELLS}

Effective AssemblyType (on the subgrid) for two AssemblyTypes where the first one is related to where the finite element functions live and the second one to where something should be assembled both with respect to the common parent grid (e.g. face-based finite elements live on a subgrid of all faces, where the faces are the cells in this subgrid, and they cannot be evaluated over the cells of the parentgrid, but on the faces of the parengrid, which are the cells in the subgrid)

source
ExtendableFEMBase.boundarydofsMethod
boundarydofs(FES; dofmap, regions) -> Any

returns an array with the number of all dofs associated to a boundary dofmap (default: BFaceDofs) and certain boundary regions (default: all regions)

source
ExtendableFEMBase.get_ndofsMethod
get_ndofs(
    P::ExtendableFEMBase.ParsedDofMap,
    _::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64

total number of dofs of the DofType

source
ExtendableFEMBase.get_ndofs4cMethod
get_ndofs4c(
    P::ExtendableFEMBase.ParsedDofMap,
    _::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64

total number of dofs of the DofType for a single component

source
ExtendableFEMBase.init_broken_dofmap!Method
init_broken_dofmap!(
    FES::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    DM::Union{Type{BFaceDofs}, Type{FaceDofs}}
) -> ExtendableGrids.VariableTargetAdjacency

generates the BFaceDofs or FaceDofs DofMap for a broken FESpace based on the pattern of the FEType

source
ExtendableFEMBase.init_dofmap_from_pattern!Method
init_dofmap_from_pattern!(
    FES::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
    DM::Type{<:DofMap}
) -> Any

generates the DofMap for the FESpace based on the pattern of the FEType

source
ExtendableFEMBase.parse_patternMethod
parse_pattern(
    pattern::String
) -> Vector{ExtendableFEMBase.DofMapPatternSegment}

splits a pattern string into pairs of single chars and Ints and generated a sequence of DofMapPatternSegments

source

The following DofMap subtypes are available and are used as keys to access the dofmap via FESpace[DofMap] (which is equivalent to FESpace.dofmaps[DofMap]).

DofMapExplanation
CellDofsdegrees of freedom for on each cell
FaceDofsdegrees of freedom for each face
EdgeDofsdegrees of freedom for each edge (in 3D)
BFaceDofsdegrees of freedom for each boundary face
BEdgeDofsdegrees of freedom for each boundary edge (in 3D)