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.AbstractFiniteElement — Type
abstract type AbstractFiniteElementroot type for finite element types
ExtendableFEMBase.AbstractH1FiniteElement — Type
abstract type AbstractH1FiniteElement <: AbstractFiniteElementroot type for H1-conforming finite element types
ExtendableFEMBase.AbstractH1FiniteElementWithCoefficients — Type
abstract type AbstractH1FiniteElementWithCoefficients <: AbstractH1FiniteElementroot type for H1-conforming finite element types with additional coefficients
ExtendableFEMBase.AbstractHcurlFiniteElement — Type
abstract type AbstractHcurlFiniteElement <: AbstractFiniteElementroot type for Hcurl-conforming finite element types with additional coefficients
ExtendableFEMBase.AbstractHdivFiniteElement — Type
abstract type AbstractHdivFiniteElement <: AbstractFiniteElementroot type for Hdiv-conforming finite element types
ExtendableFEMBase.FESpace — Type
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}
endA 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 ofxgrid).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.
ExtendableFEMBase.FESpace — Method
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
FESpacerepresents 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
brokenswitch allows the creation of a "broken" (discontinuous) finite element space, where basis functions are not required to be continuous across element boundaries. - The
regionsargument can be used to restrict the space to a subset of the grid. - If no
ATis 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).
Base.eltype — Method
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.
Base.getindex — Method
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.
Base.setindex! — Method
setindex!(FES::FESpace, v, DM::Type{<:DofMap}) -> Any
Set new dofmap
ExtendableFEMBase.assemblytype — Method
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.
ExtendableFEMBase.broken — Method
broken(FES::FESpace) -> Bool
returns true if the finite element space is broken, false if not
ExtendableFEMBase.coffset — Method
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).
ExtendableFEMBase.count_ndofs — Method
count_ndofs(
xgrid,
FEType,
broken::Bool
) -> Tuple{Int64, Int64}
counts the total number of degrees of freedom for the FEType for the whole grid
ExtendableFEMBase.dofgrid — Method
dofgrid(FES::FESpace) -> ExtendableGrids.ExtendableGrid
returns the dofgrid of the finite element space.
ExtendableFEMBase.get_AT — Method
get_AT(_::FESpace{Tv, Ti, FEType, AT}) -> Any
returns the support of the finite element space
ExtendableFEMBase.get_FEType — Method
get_FEType(_::FESpace{Tv, Ti, FEType, AT}) -> Any
returns the finite element type of the finite element space
ExtendableFEMBase.get_basis — Method
get_basis(
AT::Type{<:ExtendableGrids.AssemblyType},
FEType::Type{<:AbstractFiniteElement},
EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#get_basis##40"
returns the closure function of form
closure(refbasis, xref)that computes the evaluations of the basis functions on the reference geometry
ExtendableFEMBase.get_basissubset — Method
get_basissubset(
::Type{<:ExtendableGrids.AssemblyType},
FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#get_basissubset##5"{_A, Int64} where _A
get_basissubset(
::Type{<:ExtendableGrids.AssemblyType},
FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
::Type{<:ExtendableGrids.AbstractElementGeometry},
xgrid
) -> Union{ExtendableFEMBase.var"#20#21", ExtendableFEMBase.var"#closure#get_basissubset##6"{_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)
ExtendableFEMBase.get_coefficients — Method
get_coefficients(
::Type{<:ExtendableGrids.AssemblyType},
FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> ExtendableFEMBase.var"#closure#get_coefficients##15"
get_coefficients(
::Type{<:ExtendableGrids.AssemblyType},
FE::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT},
::Type{<:ExtendableGrids.AbstractElementGeometry},
xgrid
) -> ExtendableFEMBase.var"#closure#get_coefficients##21"
returns the coefficients for local evaluations of finite element functions ( see e.g. h1v_br.jl for a use-case)
ExtendableFEMBase.get_edim — Method
get_edim(FEType::Type{<:AbstractFiniteElement}) -> Any
returns the element dimension of the finite element
ExtendableFEMBase.get_ncomponents — Method
get_ncomponents(FES::FESpace) -> Any
returns the number of components of the FESpace (= number of components of its FEType)
ExtendableFEMBase.get_ndofs — Method
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
ExtendableFEMBase.get_ndofs_all — Method
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
ExtendableFEMBase.get_polynomialorder — Method
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)
ExtendableFEMBase.interior_dofs_offset — Method
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)
ExtendableFEMBase.isdefined — Method
isdefined(
FEType::Type{<:AbstractFiniteElement},
EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Bool
tells if FEType is defined on this ElementGeometry
ExtendableFEMBase.name — Method
name(FES::FESpace) -> String
returns the name of the finite element space.
ExtendableFEMBase.ndofs — Method
ndofs(FES::FESpace) -> Int64
returns the total number of degrees of freedom of the finite element space.
ExtendableFEMBase.xgrid — Method
xgrid(FES::FESpace) -> ExtendableGrids.ExtendableGrid
returns the computational grid of the finite element space.
DofMaps
ExtendableFEMBase.BEdgeDofs — Type
abstract type BEdgeDofs <: DofMapKey type describing the dofs for each boundary edge of the dofgrid
ExtendableFEMBase.BEdgeDofsParent — Type
abstract type BEdgeDofsParent <: DofMapKey type describing the dofs for each boundary edge of the parentgrid
ExtendableFEMBase.BFaceDofs — Type
abstract type BFaceDofs <: DofMapKey type describing the dofs for each boundary face of the dofgrid
ExtendableFEMBase.BFaceDofsParent — Type
abstract type BFaceDofsParent <: DofMapKey type describing the dofs for each boundary face of the parentgrid
ExtendableFEMBase.CellDofs — Type
abstract type CellDofs <: DofMapKey type describing the dofs for each cell of the dofgrid
ExtendableFEMBase.CellDofsParent — Type
abstract type CellDofsParent <: DofMapKey type describing the dofs for each cell of the parentgrid
ExtendableFEMBase.DofMap — Type
abstract type DofMap <: ExtendableGrids.AbstractGridAdjacencyA 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
DofMapsubtypes (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 mapExtendableFEMBase.DofMapPatternSegment — Type
struct DofMapPatternSegmentPattern segment of a dofmap sequence
ExtendableFEMBase.DofType — Type
abstract type DofTypeAbstract type for all dof types
ExtendableFEMBase.DofType — Method
DofType(c::Char) -> Union{Nothing, Type}
parses a Char into a DofType
ExtendableFEMBase.DofTypeEdge — Type
abstract type DofTypeEdge <: ExtendableFEMBase.DofTypeDof connected to an edge
ExtendableFEMBase.DofTypeFace — Type
abstract type DofTypeFace <: ExtendableFEMBase.DofTypeDof connected to a face
ExtendableFEMBase.DofTypeInterior — Type
abstract type DofTypeInterior <: ExtendableFEMBase.DofTypeDof connected to the interior of an item
ExtendableFEMBase.DofTypeNode — Type
abstract type DofTypeNode <: ExtendableFEMBase.DofTypeDof connected to a single vertex
ExtendableFEMBase.DofTypePCell — Type
abstract type DofTypePCell <: ExtendableFEMBase.DofTypeDof connected to a parent cell
ExtendableFEMBase.EdgeDofs — Type
abstract type EdgeDofs <: DofMapKey type describing the dofs for each edge of the dofgrid
ExtendableFEMBase.EdgeDofsParent — Type
abstract type EdgeDofsParent <: DofMapKey type describing the dofs for each edge of the parentgrid
ExtendableFEMBase.FaceDofs — Type
abstract type FaceDofs <: DofMapKey type describing the dofs for each face of the dofgrid
ExtendableFEMBase.FaceDofsParent — Type
abstract type FaceDofsParent <: DofMapKey type describing the dofs for each face of the parentgrid
ExtendableFEMBase.ParsedDofMap — Type
struct ParsedDofMapParsed dofmap, basically a sequence of DofMapPatternSegment
ExtendableFEMBase.ParsedDofMap — Method
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
ExtendableFEMBase.Dofmap4AssemblyType — Method
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)
ExtendableFEMBase.Dofmap4AssemblyType — Method
Dofmap4AssemblyType(_::Type{ON_CELLS}) -> Type{CellDofs}
Default Dofmap for AssemblyType
ExtendableFEMBase.EffAT4AssemblyType — Method
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)
ExtendableFEMBase.ItemEdges4DofMap — Method
ItemEdges4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellEdges}
ItemEdges grid component for dofmaps
ExtendableFEMBase.ItemGeometries4DofMap — Method
ItemGeometries4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellGeometries}
ItemGeomtries grid component for dofmaps
ExtendableFEMBase.ParentDofmap4Dofmap — Method
ParentDofmap4Dofmap(
_::Type{CellDofs}
) -> Type{CellDofsParent}
Parent Dofmap for Dofmap
ExtendableFEMBase.Sub2Sup4DofMap — Method
Sub2Sup4DofMap(
_::Type{<:DofMap}
) -> Type{ExtendableGrids.BFaceFaces}
SubItemEdges grid component for dofmaps
ExtendableFEMBase.SuperItemNodes4DofMap — Method
SuperItemNodes4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellNodes}
ItemNodes grid component for dofmaps
ExtendableFEMBase.UCG4DofMap — Method
UCG4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.UniqueCellGeometries}
Unique geometry grid component for dofmaps
ExtendableFEMBase.boundarydofs — Method
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)
ExtendableFEMBase.get_ndofs — Method
get_ndofs(
P::ExtendableFEMBase.ParsedDofMap,
_::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64
total number of dofs of the DofType
ExtendableFEMBase.get_ndofs4c — Method
get_ndofs4c(
P::ExtendableFEMBase.ParsedDofMap,
_::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64
total number of dofs of the DofType for a single component
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
ExtendableFEMBase.init_dofmap! — Method
init_dofmap!(FES::FESpace, DM::Type{<:DofMap}) -> Any
generates the requested DofMap for the FESpace
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
ExtendableFEMBase.parse_pattern — Method
parse_pattern(
pattern::String
) -> Vector{ExtendableFEMBase.DofMapPatternSegment}
splits a pattern string into pairs of single chars and Ints and generated a sequence of DofMapPatternSegments
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]).
| DofMap | Explanation |
|---|---|
| CellDofs | degrees of freedom for on each cell |
| FaceDofs | degrees of freedom for each face |
| EdgeDofs | degrees of freedom for each edge (in 3D) |
| BFaceDofs | degrees of freedom for each boundary face |
| BEdgeDofs | degrees of freedom for each boundary edge (in 3D) |