Misc methods
Test matrix generation
ExtendableSparse.fdrand! — Method
fdrand!(A, nx; ...)
fdrand!(A, nx, ny; ...)
fdrand!(A, nx, ny, nz; update, rand)
After setting all nonzero entries to zero, fill resp. update matrix with finite difference discretization data on a unit hypercube. See fdrand for documentation of the parameters.
It is required that size(A)==(N,N) where N=nx*ny*nz
ExtendableSparse.fdrand — Method
fdrand(, nx; ...)
fdrand(, nx, ny; ...)
fdrand(, nx, ny, nz; matrixtype, update, rand, symmetric)
fdrand(nx; ...)
Create matrix for a mock finite difference operator for a diffusion problem with random coefficients on a unit hypercube $\Omega\subset\mathbb R^d$. with $d=1$ if nx>1 && ny==1 && nz==1, $d=2$ if nx>1 && ny>1 && nz==1 and $d=3$ if nx>1 && ny>1 && nz>1 . In the symmetric case it corresponds to
\[ \begin{align*} -\nabla a \nabla u &= f&& \text{in}\; \Omega \\ a\nabla u\cdot \vec n + b u &=g && \text{on}\; \partial\Omega \end{align*}\]
The matrix is irreducibly diagonally dominant, has positive main diagonal entries and nonpositive off-diagonal entries, hence it has the M-Property. Therefore, its inverse will be a dense matrix with positive entries, and the spectral radius of the Jacobi iteration matrix $ho(I-D(A)^{-1}A)<1$ .
Moreover, in the symmetric case, it is positive definite.
Parameters+ default values:
| Parameter + default vale | Description |
|---|---|
nx | Number of unknowns in x direction |
ny | Number of unknowns in y direction |
nz | Number of unknowns in z direction |
matrixtype = SparseMatrixCSC | Matrix type |
update = (A,v,i,j)-> A[i,j]+=v | Element update function |
rand =()-> rand() | Random number generator |
symmetric=true | Whether to create symmetric matrix or not |
The sparsity structure is fixed to an orthogonal grid, resulting in a 3, 5 or 7 diagonal matrix depending on dimension. The entries are random unless e.g. rand=()->1 is passed as random number generator. Tested for Matrix, SparseMatrixCSC, ExtendableSparseMatrix, Tridiagonal, SparseMatrixLNK and :COO
ExtendableSparse.fdrand_coo — Function
fdrand_coo(T, nx; ...)
fdrand_coo(T, nx, ny; ...)
fdrand_coo(T, nx, ny, nz; rand)
Create SparseMatrixCSC via COO intermedite arrays. Just for benchmarking.
ExtendableSparse.sprand! — Method
sprand!(A, xnnz)
Fill empty sparse matrix A with random nonzero elements from interval [1,2] using incremental assembly.
ExtendableSparse.sprand_sdd! — Method
sprand_sdd!(A; nnzrow)
Fill sparse matrix with random entries such that it becomes strictly diagonally dominant and thus invertible and has a fixed number nnzrow (default: 4) of nonzeros in its rows. The matrix bandwidth is bounded by sqrt(n) in order to resemble a typical matrix of a 2D piecewise linear FEM discretization.
Methods for SparseMatrixCSC
SparseArrays.SparseMatrixCSC — Method
SparseArrays.SparseMatrixCSC(A::AbstractExtendableSparseMatrixCSC)Create SparseMatrixCSC from ExtendableSparseMatrix
Base.:\ — Method
Base.:\(::AbstractExtendableSparseMatrixCSC, b)\ for ExtendableSparse. It calls the LU factorization form Sparspak.jl, unless GPL components are allowed in the Julia sysimage and the floating point type of the matrix is Float64 or Complex64. In that case, Julias standard `` is called, which is realized via UMFPACK.
Base.eltype — Method
Base.eltype(::AbstractExtendableSparseMatrixCSC{Tv, Ti})Return element type.
Base.getindex — Method
getindex(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC,
i::Integer,
j::Integer
) -> Any
Base.getindex — Method
getindex(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
i::Integer,
j::Integer
) -> Any
Base.setindex! — Method
setindex!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC,
v,
i::Integer,
j::Integer
) -> Any
Base.setindex! — Method
setindex!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
v,
i::Integer,
j::Integer
) -> AbstractVecOrMat
Base.similar — Method
similar(
m::ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti}
) -> Any
Base.similar — Method
similar(
m::ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti},
_::Type{T}
) -> Any
ExtendableSparse.eliminate_dirichlet! — Method
eliminate_dirichlet!(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
dirichlet
) -> ExtendableSparse.AbstractExtendableSparseMatrixCSC
ExtendableSparse.eliminate_dirichlet! — Method
eliminate_dirichlet!(A,dirichlet_marker)Eliminate dirichlet nodes in matrix by setting
A[:,i]=0; A[i,:]=0; A[i,i]=1for a node i marked as Dirichlet.
Returns A.
ExtendableSparse.eliminate_dirichlet — Method
eliminate_dirichlet(A,dirichlet_marker)Create a copy B of A sharing the sparsity pattern. Eliminate dirichlet nodes in B by setting
B[:,i]=0; B[i,:]=0; B[i,i]=1for a node i marked as Dirichlet.
Returns B.
ExtendableSparse.eliminate_dirichlet — Method
eliminate_dirichlet(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
dirichlet
) -> ExtendableSparse.GenericExtendableSparseMatrixCSC
ExtendableSparse.findindex — Method
findindex(csc::SparseMatrixCSC{T}, i, j) -> Int64
Return index corresponding to entry [i,j] in the array of nonzeros, if the entry exists, otherwise, return 0.
ExtendableSparse.flush! — Method
flush!(csc::SparseMatrixCSC) -> SparseMatrixCSC
Trivial flush! method for allowing to run the code with both ExtendableSparseMatrix and SparseMatrixCSC.
ExtendableSparse.flush! — Method
flush!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti}
) -> ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti} where {Tm, Tv, Ti}
ExtendableSparse.flush! — Method
flush!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti}
) -> ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti} where {Tm, Tv, Ti}
ExtendableSparse.mark_dirichlet — Method
mark_dirichlet(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC;
penalty
) -> Vector{Bool}
ExtendableSparse.mark_dirichlet — Method
mark_dirichlet(A; penalty=1.0e20)
Return boolean vector marking Dirichlet nodes, known by A[i,i]>=penalty
ExtendableSparse.nnznew — Method
nnznew(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC
) -> Integer
ExtendableSparse.nnznew — Method
nnznew(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC
) -> Real
ExtendableSparse.partitioning! — Method
partitioning!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti},
colparts,
partnodes
) -> ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti} where {Tm, Tv, Ti}
Set node partitioning.
ExtendableSparse.pattern_equal — Method
pattern_equal(a::SparseMatrixCSC,b::SparseMatrixCSC)Check if sparsity patterns of two SparseMatrixCSC objects are equal. This is generally faster than comparing hashes.
ExtendableSparse.phash — Method
phash(csc::SparseMatrixCSC) -> Any
Hash of csc matrix pattern.
ExtendableSparse.rawupdateindex! — Function
rawupdateindex!(A::AbstractExtendableSparseMatrixCSC,op,v,i,j,part = 1)Add or update entry of A: A[i,j]=op(A[i,j],v) without checking if a zero is inserted. The optional parameter part denotes the partition.
ExtendableSparse.rawupdateindex! — Function
rawupdateindex!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
op,
v,
i,
j
) -> Any
rawupdateindex!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
op,
v,
i,
j,
tid
) -> Any
ExtendableSparse.rawupdateindex! — Function
rawupdateindex!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC,
op,
v,
i,
j
) -> Any
rawupdateindex!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC,
op,
v,
i,
j,
part
) -> Any
ExtendableSparse.reset! — Method
reset!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC
) -> ExtendableSparse.GenericMTExtendableSparseMatrixCSC
ExtendableSparse.reset! — Method
reset!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti}
) -> ExtendableSparse.GenericExtendableSparseMatrixCSC{Tm, Tv, Ti} where {Tm, Tv, Ti}
ExtendableSparse.reset! — Method
reset!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti},
p::Integer
) -> ExtendableSparse.GenericMTExtendableSparseMatrixCSC{Tm, Tv, Ti} where {Tm, Tv, Ti}
ExtendableSparse.updateindex! — Function
updateindex!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
op,
v,
i,
j
) -> Any
updateindex!(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
op,
v,
i,
j,
tid
) -> Any
ExtendableSparse.updateindex! — Method
updateindex!(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC,
op,
v,
i,
j
) -> Any
ExtendableSparse.updateindex! — Method
updateindex!(
csc::SparseMatrixCSC{Tv, Ti<:Integer},
op,
v,
i,
j
) -> SparseMatrixCSC
Update element of the matrix with operation op without introducing new nonzero elements.
This can replace the following code and save one index search during access:
using ExtendableSparse # hide
using SparseArrays # hide
A=spzeros(3,3)
A[1,2]+=0.1
Ausing ExtendableSparse # hide
using SparseArrays # hide
A=spzeros(3,3)
updateindex!(A,+,0.1,1,2)
ALinearAlgebra.cond — Function
cond(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC
) -> Any
cond(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
p::Real
) -> Any
flush! and calculate cond from cscmatrix
LinearAlgebra.issymmetric — Method
issymmetric(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC
) -> Bool
flush! and check for symmetry of cscmatrix
LinearAlgebra.ldiv! — Method
ldiv!(
r::AbstractArray,
ext::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
x::AbstractArray
) -> SparseMatrixCSC
flush! and ldiv! with ext.cscmatrix
LinearAlgebra.mul! — Method
mul!(
r::AbstractVecOrMat,
ext::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
x::AbstractVecOrMat
)
flush! and multiply with ext.cscmatrix
LinearAlgebra.mul! — Method
mul!(
r::AbstractVecOrMat,
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC,
x::AbstractVecOrMat
)
LinearAlgebra.norm — Function
norm(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC
) -> Any
norm(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
p::Real
) -> Any
flush! and calculate norm from cscmatrix
LinearAlgebra.opnorm — Function
opnorm(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC
) -> Any
opnorm(
A::ExtendableSparse.AbstractExtendableSparseMatrixCSC,
p::Real
) -> Any
flush! and calculate opnorm from cscmatrix
SparseArrays.dropzeros! — Method
dropzeros!(
ext::ExtendableSparse.AbstractExtendableSparseMatrixCSC
) -> SparseMatrixCSC
SparseArrays.findnz — Method
SparseArrays.findnz(ext::AbstractExtendableSparseMatrixCSC)flush! and return findnz(ext.cscmatrix).
SparseArrays.getcolptr — Method
SparseArrays.getcolptr(ext::AbstractExtendableSparseMatrixCSC)flush! and return colptr of in ext.cscmatrix.
SparseArrays.nnz — Method
SparseArrays.nnz(ext::AbstractExtendableSparseMatrixCSC)flush! and return number of nonzeros in ext.cscmatrix.
SparseArrays.nonzeros — Method
SparseArrays.nonzeros(ext::AbstractExtendableSparseMatrixCSC) = nonzeros(sparse(ext))flush! and return nonzeros in ext.cscmatrix.
SparseArrays.rowvals — Method
SparseArrays.rowvals(ext::AbstractExtendableSparseMatrixCSC)flush! and return rowvals in ext.cscmatrix.
SparseArrays.sparse — Method
SparseArrays.sparse(A::AbstractExtendableSparseMatrixCSC)Return SparseMatrixCSC which contains all matrix entries introduced so far.
SparseArrays.sparse — Method
sparse(
ext::ExtendableSparse.GenericExtendableSparseMatrixCSC
) -> SparseMatrixCSC
SparseArrays.sparse — Method
sparse(
ext::ExtendableSparse.GenericMTExtendableSparseMatrixCSC
) -> SparseMatrixCSC