FESpace

To generate a finite element space only a finite element type and a grid is needed, dofmaps are generated automatically on their first demand.

ExtendableFEMBase.FESpaceType
struct FESpace{Tv, Ti, FEType<:AbstractFiniteElement,AT<:AssemblyType}
	name::String                          # full name of finite element space (used in messages)
	broken::Bool                          # if true, broken dofmaps are generated
	ndofs::Int                            # total number of dofs
	coffset::Int                          # offset for component dofs
	xgrid::ExtendableGrid[Tv,Ti}          # link to xgrid 
	dofgrid::ExtendableGrid{Tv,Ti}	      # link to (sub) grid used for dof numbering (expected to be equal to or child grid of xgrid)
	dofmaps::Dict{Type{<:AbstractGridComponent},Any} # backpack with dofmaps
end

A struct that has a finite element type as parameter and carries dofmaps (CellDofs, FaceDofs, BFaceDofs) plus additional grid information and access to arrays holding coefficients if needed.

source
ExtendableFEMBase.FESpaceMethod
function FESpace{FEType<:AbstractFiniteElement,AT<:AssemblyType}(
	xgrid::ExtendableGrid{Tv,Ti};
	name = "",
	broken::Bool = false)

Constructor for FESpace of the given FEType, AT = ONCELLS/ONFACES/ONEDGES generates a finite elements space on the cells/faces/edges of the provided xgrid (if omitted ONCELLS is used as default). The broken switch allows to generate a broken finite element space (that is piecewise H1/Hdiv/HCurl). If no name is provided it is generated automatically from FEType. If no AT is provided, the space is generated ON_CELLS.

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.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#123"

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}
) -> Union{ExtendableFEMBase.var"#21#22", ExtendableFEMBase.var"#closure#268"{_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#195"

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}
) -> Int64

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

DofMaps

ExtendableFEMBase.DofMapType
abstract type DofMap <: ExtendableGrids.AbstractGridAdjacency

Dofmaps are stored as an ExtendableGrids.AbstractGridAdjacency in the finite element space and collect information with respect to different AssemblyTypes. They are generated automatically on demand and the dofmaps associated to each subtype can be accessed via FESpace[DofMap].

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)