FEMatrix
A FEMatrix
represents a block-structured finite element matrix, where each block (a FEMatrixBlock
) corresponds to a pair of finite element spaces and operates on a submatrix of a shared ExtendableSparseMatrix
.
ExtendableFEMBase.FEMatrix
— TypeFEMatrix(FESX::Union{FESpace, Vector{<:FESpace}}, FESY::Union{FESpace, Vector{<:FESpace}}=FESX; TvM=Float64, TiM=Int64, ...)
Constructs an FEMatrix
for storing (sparse) matrix representations associated with one or more pairs of finite element spaces (FESpace
).
– If FESX
and FESY
are single FESpace
objects, the resulting FEMatrix
contains one rectangular block.
- If
FESX
and/orFESY
are vectors ofFESpace
objects, the resultingFEMatrix
is block-structured, with one block for each pair.
Optionally, you can assign a name (as a String
for all blocks) and/or tags (as arrays for rows and columns) to the blocks for identification and access.
Arguments
FESX::FESpace
orFESX::Vector{<:FESpace}
: Row finite element space(s).FESY::FESpace
orFESY::Vector{<:FESpace}
: Column finite element space(s) (default: same as FESX).
Keyword Arguments
entries
: Optional sparse matrix of coefficients. If not provided, a new sparse matrix of appropriate size is created.name
: Name for the matrix or for each block (default::automatic
).tags
: Array of tags for both rows and columns (default:nothing
).tagsX
: Array of tags for the row blocks (default:tags
).tagsY
: Array of tags for the column blocks (default:tagsX
).npartitions
: Number of partitions for the underlying sparse matrix (default:1
).- Additional keyword arguments are passed to the underlying block constructors.
Returns
- An
FEMatrix
object with one or moreFEMatrixBlock
s, each corresponding to a given pair ofFESpace
objects.
ExtendableFEMBase.FEMatrix
— Typestruct FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal} <: SparseArrays.AbstractSparseArray{TvM, TiM, 2}
A block-structured sparse matrix type for finite element assembly.
FEMatrix
represents a (typically sparse) matrix with an additional layer of structure: it is subdivided into multiple FEMatrixBlock
s, each associated with a specific pair of finite element spaces (FESpace
).
Fields
FEMatrixBlocks::Array{FEMatrixBlock{TvM, TiM, TvG, TiG}, 1}
: The array of matrix blocks, each corresponding to a pair of row and column finite element spaces.entries::AbstractSparseMatrix{TvM, TiM}
: The underlying sparse matrix storing all coefficients.tags::Matrix{Any}
: Optional tags for identifying or accessing blocks.
Type Parameters
TvM
: Value type for matrix entries (e.g.,Float64
).TiM
: Integer type for matrix indices (e.g.,Int64
).TvG
,TiG
: Value and index types for the associated finite element spaces.nbrow
: Number of block rows.nbcol
: Number of block columns.nbtotal
: Total number of blocks.
ExtendableFEMBase.FEMatrixBlock
— Typestruct FEMatrixBlock{TvM, TiM, TvG, TiG, FETypeX, FETypeY, APTX, APTY} <: AbstractArray{TvM, 2}
A block of an FEMatrix
representing the coupling between two finite element spaces.
FEMatrixBlock
acts as a two-dimensional array (subclassing AbstractArray{TvM, 2}
) and stores the coefficients for a specific pair of row and column finite element spaces (FESpace
). Each block is mapped to a submatrix of the global sparse matrix, with offsets and sizes corresponding to the degrees of freedom of the associated spaces.
Fields
name::String
: Name of the block (for identification and display).FES::FESpace{TvG, TiG, FETypeX, APTX}
: Row finite element space.FESY::FESpace{TvG, TiG, FETypeY, APTY}
: Column finite element space.offset::Int64
: Row offset in the global matrix.offsetY::Int64
: Column offset in the global matrix.last_index::Int64
: Last row index for this block.last_indexY::Int64
: Last column index for this block.entries::AbstractSparseMatrix{TvM, TiM}
: Reference to the underlying global sparse matrix (shared with the parentFEMatrix
).
See also: [FEMatrix
], [FESpace
]
Base.fill!
— Methodfill!(B::FEMatrixBlock{Tv, Ti}, value)
Custom fill
function for FEMatrixBlock
(only fills the already present nzval in 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, FEB::FEMatrixBlock)
Custom show
function for FEMatrixBlock
that prints its coordinates and the name.
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
ExtendableFEMBase._addnz
— Methodadds value v to matrix at position (i,j) if it is nonzero
ExtendableFEMBase.add!
— Methodadd!(A::FEMatrix{Tv, Ti}, B::FEMatrix{Tv, Ti}; kwargs...)
Adds FEMatrix/ExtendableSparseMatrix/CSCMatrix B to FEMatrix A.
ExtendableFEMBase.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv, Ti},
B::FEMatrixBlock{Tv, Ti};
factor,
transpose
)
Adds FEMatrixBlock B to FEMatrixBlock A.
ExtendableFEMBase.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer};
factor,
transpose
)
Adds ExtendableSparseMatrix B to FEMatrixBlock A.
ExtendableFEMBase.addblock!
— Methodaddblock!(
A::FEMatrixBlock{Tv},
cscmat::SparseArrays.SparseMatrixCSC{Tv, Ti<:Integer};
factor,
transpose
)
Adds SparseMatrixCSC B to FEMatrixBlock A.
ExtendableFEMBase.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.
ExtendableFEMBase.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.
ExtendableFEMBase.addblock_matmul!
— Methodaddblock_matmul!(
a::FEVectorBlock{Tv},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer},
b::FEVectorBlock{Tv};
factor
)
Adds matrix-vector product B times b to FEVectorBlock a.
ExtendableFEMBase.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.
ExtendableFEMBase.apply_nonzero_pattern!
— Methodpre-allocates the expected pattern for the default dofmaps for the AssemblyType
ExtendableFEMBase.apply_penalties!
— Methodapply_penalties!(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
fixed_dofs,
penalty
)
sets penalty to the diagonal entries of fixed_dofs in A
ExtendableFEMBase.ldrdmatmul
— Methodldrdmatmul(
a1::AbstractArray{Tv, 1},
a2::AbstractArray{Tv, 1},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer},
b1::AbstractArray{Tv, 1},
b2::AbstractArray{Tv, 1};
factor
) -> Any
Computes vector'-matrix-vector product (a1-a2)'B(b1-b2).
ExtendableFEMBase.lrmatmul
— Methodlrmatmul(
a::AbstractArray{Tv, 1},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer},
b::AbstractArray{Tv, 1};
factor
) -> Any
Computes vector'-matrix-vector product a'Bb.
ExtendableFEMBase.nbcols
— Methodnbcols(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each row.
ExtendableFEMBase.nbrows
— Methodnbrows(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each column.
ExtendableFEMBase.submatrix
— Methodsubmatrix(
A::FEMatrixBlock{Tv, Ti}
) -> ExtendableSparse.ExtendableSparseMatrixCSC
Returns the FEMatrixBlock as an ExtendableSparseMatrix
ExtendableFEMBase.submatrix
— Methodsubmatrix(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti},
srows,
scols;
factor
) -> ExtendableSparse.ExtendableSparseMatrixCSC
Generates an ExtendableSparseMatrix from the submatrix for the given row and col numbers