Finite Element Spaces and Arrays
This page describes the structure FESpace that acts as a finite element space on a given grid and provides the associated degree of freedom maps DofMaps on demand. See Implemented Finite Elements for a list of available finite element types.
Moreover, there are special arrays FEVector and FEMatrix that carry coefficients and discretised PDEOperators.
FESpace
To generate a finite element space only a finite element type and a grid is needed, dofmaps are generated automatically on demand.
GradientRobustMultiPhysics.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
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.
GradientRobustMultiPhysics.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.
GradientRobustMultiPhysics.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.
DofMaps
GradientRobustMultiPhysics.DofMap
— Typeabstract type DofMap <: 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].
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) |
FEVector
A FEVector consists of FEVectorBlocks that share a common one-dimensional array. Each block is associated to a FESpace and can only write into a region of the common array specified by offsets. It also acts as a one-dimensional AbstractArray itself.
GradientRobustMultiPhysics.FEVector
— Typestruct FEVector{T, Tv, Ti}
a plain array but with an additional layer of several FEVectorBlock subdivisions each carrying coefficients for their associated FESpace
GradientRobustMultiPhysics.FEVector
— MethodFEVector{T}(FES; name = "auto") where T <: Real
Creates FEVector that has one block if FES is a single FESpace, and a blockwise FEVector if FES is a vector of FESpaces. Optionally a name for the vector (as a String), or each of the blocks (as a vector of Strings) can be specified.
GradientRobustMultiPhysics.FEVectorBlock
— Typestruct FEVectorBlock{T, Tv, Ti, FEType, APT} <: AbstractArray{T, 1}
block of an FEVector that carries coefficients for an associated FESpace and can be assigned as an AbstractArray (getindex, setindex, size, length)
Base.append!
— Methodappend!(
FEF::FEVector{T},
name::String,
FES::FESpace{Tv, Ti, FEType, APT}
)
Custom append
function for FEVector
that adds a FEVectorBlock at the end.
Base.fill!
— Methodfill!(b::FEVectorBlock, value)
Custom fill
function for FEVectorBlock
(only fills the block, not the complete FEVector).
Base.length
— Methodlength(FEB::FEVectorBlock) -> Int64
Custom length
function for FEVectorBlock
that gives the coressponding number of degrees of freedoms of the associated FESpace
Base.length
— Methodlength(FEF::FEVector) -> Int64
Custom length
function for FEVector
that gives the number of defined FEMatrixBlocks in it
Base.show
— Methodshow(io::IO, FEF::FEVector)
Custom show
function for FEVector
that prints some information on its blocks.
GradientRobustMultiPhysics.FESpaces
— MethodFESpaces(FEV::FEVector{T, Tv, Ti}) -> Any
Returns the vector of FEspaces for the blocks of the given FEVector.
GradientRobustMultiPhysics.addblock!
— Methodaddblock!(a::FEVectorBlock, b::AbstractVector; factor)
Adds Array b to FEVectorBlock a.
GradientRobustMultiPhysics.addblock!
— Methodaddblock!(a::FEVectorBlock, b::FEVectorBlock; factor)
Adds FEVectorBlock b to FEVectorBlock a.
LinearAlgebra.dot
— Methoddot(a::FEVectorBlock{T}, b::FEVectorBlock{T}) -> Any
Scalar product between two FEVEctorBlocks
FEMatrix
A FEMatrix consists of FEMatrixBlocks that share a common ExtendableSparseMatrix. Each block is associated to two FESpaces and can only write into a submatrix of the common sparse matrix specified by offsets. It also acts as a two-dimensional AbstractArray itself.
GradientRobustMultiPhysics.FEMatrix
— Typestruct FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
an AbstractMatrix (e.g. an ExtendableSparseMatrix) with an additional layer of several FEMatrixBlock subdivisions each carrying coefficients for their associated pair of FESpaces
GradientRobustMultiPhysics.FEMatrix
— MethodFEMatrix{TvM,TiM}(FESX, FESY; name = "auto")
Creates FEMatrix with one rectangular block (FESX,FESY) if FESX and FESY are single FESpaces, or a rectangular block matrix with blocks corresponding to the entries of the FESpace vectors FESX and FESY. Optionally a name for the matrix can be given.
GradientRobustMultiPhysics.FEMatrix
— MethodFEMatrix{TvM,TiM}(name::String, FES::FESpace{TvG,TiG,FETypeX,APTX}) where {TvG,TiG,FETypeX,APTX}
Creates FEMatrix with one square block (FES,FES).
GradientRobustMultiPhysics.FEMatrixBlock
— Typestruct FEMatrixBlock{TvM, TiM, TvG, TiG, FETypeX, FETypeY, APTX, APTY} <: AbstractArray{TvM, 2}
block of an FEMatrix that carries coefficients for an associated pair of FESpaces and can be assigned as an two-dimensional AbstractArray (getindex, setindex, size)
Base.fill!
— Methodfill!(B::FEMatrixBlock{Tv, Ti}, value)
Custom fill
function for FEMatrixBlock
(only fills the block, not the complete FEMatrix).
Base.length
— Methodlength(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Custom length
function for FEMatrix
that gives the total number of defined FEMatrixBlocks in it
Base.show
— Methodshow(
io::IO,
FEM::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
)
Custom show
function for FEMatrix
that prints some information on its blocks.
Base.size
— Methodsize(FEB::FEMatrixBlock) -> Tuple{Int64, Int64}
Custom size
function for FEMatrixBlock
that gives a tuple with the size of the block (that coressponds to the number of degrees of freedoms in X and Y)
Base.size
— Methodsize(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Tuple{Any, Any}
Custom size
function for FEMatrix
that gives a tuple with the number of rows and columns of the FEBlock overlay
GradientRobustMultiPhysics.add!
— Methodadd!(
A::FEMatrix{Tv, Ti},
B::FEMatrix{Tv, Ti};
factor,
transpose
)
Adds FEMatrix B to FEMatrix A.
GradientRobustMultiPhysics.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv, Ti},
B::FEMatrixBlock{Tv, Ti};
factor,
transpose
)
Adds FEMatrixBlock B to FEMatrixBlock A.
GradientRobustMultiPhysics.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv},
B::ExtendableSparseMatrix{Tv, Ti<:Integer};
factor,
transpose
)
Adds ExtendableSparseMatrix B to FEMatrixBlock A.
GradientRobustMultiPhysics.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv},
cscmat::SparseArrays.SparseMatrixCSC{Tv, Ti<:Integer};
factor,
transpose
)
Adds SparseMatrixCSC B to FEMatrixBlock A.
GradientRobustMultiPhysics.addblock_matmul!
— Methodaddblock_matmul!(
a::AbstractArray{Tv, 1},
B::FEMatrixBlock{Tv, Ti},
b::AbstractArray{Tv, 1};
factor,
transposed
)
Adds matrix-vector product B times b to FEVectorBlock a.
GradientRobustMultiPhysics.addblock_matmul!
— Methodaddblock_matmul!(
A::FEMatrixBlock{Tv},
cscmatB::SparseArrays.SparseMatrixCSC{Tv, Ti},
cscmatC::SparseArrays.SparseMatrixCSC{Tv, Ti};
factor,
transposed
)
Adds matrix-matrix product B times C to FEMatrixBlock A.
GradientRobustMultiPhysics.addblock_matmul!
— Methodaddblock_matmul!(
a::FEVectorBlock{Tv},
B::ExtendableSparseMatrix{Tv, Ti<:Integer},
b::FEVectorBlock{Tv};
factor
)
Adds matrix-vector product B times b to FEVectorBlock a.
GradientRobustMultiPhysics.addblock_matmul!
— Methodaddblock_matmul!(
a::FEVectorBlock{Tv},
B::FEMatrixBlock{Tv, Ti},
b::FEVectorBlock{Tv};
factor,
transposed
)
Adds matrix-vector product B times b (or B' times b if transposed = true) to FEVectorBlock a.
GradientRobustMultiPhysics.ldrdmatmul
— Methodldrdmatmul(
a1::AbstractArray{Tv, 1},
a2::AbstractArray{Tv, 1},
B::ExtendableSparseMatrix{Tv, Ti<:Integer},
b1::AbstractArray{Tv, 1},
b2::AbstractArray{Tv, 1};
factor
) -> Any
Computes vector'-matrix-vector product (a1-a2)'B(b1-b2).
GradientRobustMultiPhysics.lrmatmul
— Methodlrmatmul(
a::AbstractArray{Tv, 1},
B::ExtendableSparseMatrix{Tv, Ti<:Integer},
b::AbstractArray{Tv, 1};
factor
) -> Any
Computes vector'-matrix-vector product a'Bb.
GradientRobustMultiPhysics.nbcols
— Methodnbcols(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each row.
GradientRobustMultiPhysics.nbrows
— Methodnbrows(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each column.