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 — Type
FEMatrix(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
FESXand/orFESYare vectors ofFESpaceobjects, the resultingFEMatrixis 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::FESpaceorFESX::Vector{<:FESpace}: Row finite element space(s).FESY::FESpaceorFESY::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
FEMatrixobject with one or moreFEMatrixBlocks, each corresponding to a given pair ofFESpaceobjects.
ExtendableFEMBase.FEMatrix — Type
struct 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 FEMatrixBlocks, 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 — Type
struct 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! — Method
fill!(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 — Method
length(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Custom length function for FEMatrix that gives the total number of defined FEMatrixBlocks in it
ExtendableFEMBase._addnz — Method
adds value v to matrix at position (i,j) if it is nonzero
sourceExtendableFEMBase.add! — Method
add!(A::FEMatrix{Tv, Ti}, B::FEMatrix{Tv, Ti}; kwargs...)
Adds FEMatrix/ExtendableSparseMatrix/CSCMatrix B to FEMatrix A.
sourceExtendableFEMBase.addblock! — Method
addblock!(
A::FEMatrixBlock{Tv, Ti},
B::FEMatrixBlock{Tv, Ti};
factor,
transpose
)
Adds FEMatrixBlock B to FEMatrixBlock A.
sourceExtendableFEMBase.addblock! — Method
addblock!(
A::FEMatrixBlock{Tv},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer};
factor,
transpose
)
Adds ExtendableSparseMatrix B to FEMatrixBlock A.
sourceExtendableFEMBase.addblock! — Method
addblock!(
A::FEMatrixBlock{Tv},
cscmat::SparseArrays.SparseMatrixCSC{Tv, Ti<:Integer};
factor,
transpose
)
Adds SparseMatrixCSC B to FEMatrixBlock A.
sourceExtendableFEMBase.addblock_matmul! — Method
addblock_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.
sourceExtendableFEMBase.addblock_matmul! — Method
addblock_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.
sourceExtendableFEMBase.addblock_matmul! — Method
addblock_matmul!(
a::FEVectorBlock{Tv, Tv, Ti, TVector} where {Tv, Ti, TVector<:AbstractArray{Tv, 1}},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer},
b::FEVectorBlock{Tv, Tv, Ti, TVector} where {Tv, Ti, TVector<:AbstractArray{Tv, 1}};
factor
)
Adds matrix-vector product B times b to FEVectorBlock a.
sourceExtendableFEMBase.addblock_matmul! — Method
addblock_matmul!(
a::FEVectorBlock{Tv, Tv, Ti, TVector} where {Tv, Ti, TVector<:AbstractArray{Tv, 1}},
B::FEMatrixBlock{Tv, Ti},
b::FEVectorBlock{Tv, Tv, Ti, TVector} where {Tv, Ti, TVector<:AbstractArray{Tv, 1}};
factor,
transposed
)
Adds matrix-vector product B times b (or B' times b if transposed = true) to FEVectorBlock a.
sourceExtendableFEMBase.apply_nonzero_pattern! — Method
pre-allocates the expected pattern for the default dofmaps for the AssemblyType
sourceExtendableFEMBase.apply_penalties! — Method
apply_penalties!(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
fixed_dofs,
penalty
)
sets penalty to the diagonal entries of fixed_dofs in A
sourceExtendableFEMBase.ldrdmatmul — Method
ldrdmatmul(
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).
sourceExtendableFEMBase.lrmatmul — Method
lrmatmul(
a::AbstractArray{Tv, 1},
B::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti<:Integer},
b::AbstractArray{Tv, 1};
factor
) -> Any
Computes vector'-matrix-vector product a'Bb.
sourceExtendableFEMBase.nbcols — Method
nbcols(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each row.
sourceExtendableFEMBase.nbrows — Method
nbrows(
_::FEMatrix{TvM, TiM, TvG, TiG, nbrow, nbcol, nbtotal}
) -> Any
Gives the number of FEMatrixBlocks in each column.
sourceExtendableFEMBase.submatrix — Method
submatrix(
A::FEMatrixBlock{Tv, Ti}
) -> ExtendableSparse.ExtendableSparseMatrixCSC
Returns the FEMatrixBlock as an ExtendableSparseMatrix
sourceExtendableFEMBase.submatrix — Method
submatrix(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC{Tv, Ti},
srows,
scols;
factor
) -> ExtendableSparse.ExtendableSparseMatrixCSC
Generates an ExtendableSparseMatrix from the submatrix for the given row and col numbers
source