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.AbstractFiniteElement
— Typeabstract type AbstractFiniteElement
root type for finite element types
ExtendableFEMBase.AbstractH1FiniteElement
— Typeabstract type AbstractH1FiniteElement <: AbstractFiniteElement
root type for H1-conforming finite element types
ExtendableFEMBase.AbstractH1FiniteElementWithCoefficients
— Typeabstract type AbstractH1FiniteElementWithCoefficients <: AbstractH1FiniteElement
root type for H1-conforming finite element types with additional coefficients
ExtendableFEMBase.AbstractHcurlFiniteElement
— Typeabstract type AbstractHcurlFiniteElement <: AbstractFiniteElement
root type for Hcurl-conforming finite element types with additional coefficients
ExtendableFEMBase.AbstractHdivFiniteElement
— Typeabstract type AbstractHdivFiniteElement <: AbstractFiniteElement
root type for Hdiv-conforming finite element types
ExtendableFEMBase.FESpace
— Typestruct 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.
ExtendableFEMBase.FESpace
— Methodfunction 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.
Base.eltype
— Methodeltype(
_::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.get!
— Methodget!(FES::FESpace, DM::Type{<:DofMap}) -> Any
To be called by getindex. This triggers lazy creation of non-existing dofmaps
Base.getindex
— MethodBase.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!
— Methodsetindex!(FES::FESpace, v, DM::Type{<:DofMap}) -> Any
Set new dofmap
Base.show
— Methodshow(
io::IO,
FES::FESpace{Tv, Ti, FEType<:AbstractFiniteElement, APT}
)
Custom show
function for FESpace
that prints some information and all available dofmaps.
ExtendableFEMBase.assemblytype
— Methodassemblytype(
_::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
— Methodbroken(FES::FESpace) -> Bool
returns true if the finite element space is broken, false if not
ExtendableFEMBase.count_ndofs
— Methodcount_ndofs(
xgrid,
FEType,
broken::Bool
) -> Tuple{Int64, Int64}
counts the total number of degrees of freedom for the FEType for the whole grid
ExtendableFEMBase.get_AT
— Methodget_AT(_::FESpace{Tv, Ti, FEType, AT}) -> Any
returns the support of the finite element space
ExtendableFEMBase.get_FEType
— Methodget_FEType(_::FESpace{Tv, Ti, FEType, AT}) -> Any
returns the finite element type of the finite element space
ExtendableFEMBase.get_basis
— Methodget_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
ExtendableFEMBase.get_basissubset
— Methodget_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)
ExtendableFEMBase.get_coefficients
— Methodget_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)
ExtendableFEMBase.get_edim
— Methodget_edim(FEType::Type{<:AbstractFiniteElement}) -> Any
returns the element dimension of the finite element
ExtendableFEMBase.get_ncomponents
— Methodget_ncomponents(FES::FESpace) -> Any
returns the number of components of the FESpace (= number of components of its FEType)
ExtendableFEMBase.get_ndofs
— Methodget_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
ExtendableFEMBase.get_ndofs_all
— Methodget_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
— Methodget_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
— Methodinterior_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
— Methodisdefined(
FEType::Type{<:AbstractFiniteElement},
EG::Type{<:ExtendableGrids.AbstractElementGeometry}
) -> Bool
tells if FEType is defined on this ElementGeometry
ExtendableFEMBase.ndofs
— Methodndofs(FES::FESpace) -> Int64
returns the total number of degrees of freedom of the finite element space.
DofMaps
ExtendableFEMBase.BEdgeDofs
— Typeabstract type BEdgeDofs <: DofMap
Key type describing the dofs for each boundary edge of the dofgrid
ExtendableFEMBase.BEdgeDofsParent
— Typeabstract type BEdgeDofsParent <: DofMap
Key type describing the dofs for each boundary edge of the parentgrid
ExtendableFEMBase.BFaceDofs
— Typeabstract type BFaceDofs <: DofMap
Key type describing the dofs for each boundary face of the dofgrid
ExtendableFEMBase.BFaceDofsParent
— Typeabstract type BFaceDofsParent <: DofMap
Key type describing the dofs for each boundary face of the parentgrid
ExtendableFEMBase.CellDofs
— Typeabstract type CellDofs <: DofMap
Key type describing the dofs for each cell of the dofgrid
ExtendableFEMBase.CellDofsParent
— Typeabstract type CellDofsParent <: DofMap
Key type describing the dofs for each cell of the parentgrid
ExtendableFEMBase.DofMap
— Typeabstract 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].
ExtendableFEMBase.DofMapPatternSegment
— Typestruct DofMapPatternSegment
Pattern segment of a dofmap sequence
ExtendableFEMBase.DofType
— Typeabstract type DofType
Abstract type for all dof types
ExtendableFEMBase.DofType
— MethodDofType(c::Char) -> Union{Nothing, Type}
parses a Char into a DofType
ExtendableFEMBase.DofTypeEdge
— Typeabstract type DofTypeEdge <: ExtendableFEMBase.DofType
Dof connected to an edge
ExtendableFEMBase.DofTypeFace
— Typeabstract type DofTypeFace <: ExtendableFEMBase.DofType
Dof connected to a face
ExtendableFEMBase.DofTypeInterior
— Typeabstract type DofTypeInterior <: ExtendableFEMBase.DofType
Dof connected to the interior of an item
ExtendableFEMBase.DofTypeNode
— Typeabstract type DofTypeNode <: ExtendableFEMBase.DofType
Dof connected to a single vertex
ExtendableFEMBase.DofTypePCell
— Typeabstract type DofTypePCell <: ExtendableFEMBase.DofType
Dof connected to a parent cell
ExtendableFEMBase.EdgeDofs
— Typeabstract type EdgeDofs <: DofMap
Key type describing the dofs for each edge of the dofgrid
ExtendableFEMBase.EdgeDofsParent
— Typeabstract type EdgeDofsParent <: DofMap
Key type describing the dofs for each edge of the parentgrid
ExtendableFEMBase.FaceDofs
— Typeabstract type FaceDofs <: DofMap
Key type describing the dofs for each face of the dofgrid
ExtendableFEMBase.FaceDofsParent
— Typeabstract type FaceDofsParent <: DofMap
Key type describing the dofs for each face of the parentgrid
ExtendableFEMBase.ParsedDofMap
— Typestruct ParsedDofMap
Parsed dofmap, basically a sequence of DofMapPatternSegment
ExtendableFEMBase.ParsedDofMap
— MethodParsedDofMap(
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
— MethodDofmap4AssemblyType(
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
— MethodDofmap4AssemblyType(_::Type{ON_CELLS}) -> Type{CellDofs}
Default Dofmap for AssemblyType
ExtendableFEMBase.EffAT4AssemblyType
— MethodEffAT4AssemblyType(
_::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
— MethodItemEdges4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellEdges}
ItemEdges grid component for dofmaps
ExtendableFEMBase.ItemGeometries4DofMap
— MethodItemGeometries4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellGeometries}
ItemGeomtries grid component for dofmaps
ExtendableFEMBase.ParentDofmap4Dofmap
— MethodParentDofmap4Dofmap(
_::Type{CellDofs}
) -> Type{CellDofsParent}
Parent Dofmap for Dofmap
ExtendableFEMBase.Sub2Sup4DofMap
— MethodSub2Sup4DofMap(
_::Type{<:DofMap}
) -> Type{ExtendableGrids.BEdgeEdges}
SubItemEdges grid component for dofmaps
ExtendableFEMBase.SuperItemNodes4DofMap
— MethodSuperItemNodes4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.CellNodes}
ItemNodes grid component for dofmaps
ExtendableFEMBase.UCG4DofMap
— MethodUCG4DofMap(
_::Type{CellDofs}
) -> Type{ExtendableGrids.UniqueCellGeometries}
Unique geometry grid component for dofmaps
ExtendableFEMBase.boundarydofs
— Methodboundarydofs(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
— Methodget_ndofs(
P::ExtendableFEMBase.ParsedDofMap,
_::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64
total number of dofs of the DofType
ExtendableFEMBase.get_ndofs4c
— Methodget_ndofs4c(
P::ExtendableFEMBase.ParsedDofMap,
_::Type{ExtendableFEMBase.DofTypeNode}
) -> Int64
total number of dofs of the DofType for a single component
ExtendableFEMBase.init_broken_dofmap!
— Methodinit_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!
— Methodinit_dofmap!(FES::FESpace, DM::Type{<:DofMap}) -> Any
generates the requested DofMap for the FESpace
ExtendableFEMBase.init_dofmap_from_pattern!
— Methodinit_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
— Methodparse_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) |