Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

SparseCSCInterface sparspaklu! check sparsity pattern #27

Closed
wants to merge 6 commits into from
Closed
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
88 changes: 57 additions & 31 deletions src/SparseCSCInterface/SparseCSCInterface.jl
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
module SparseCSCInterface
using SparseArrays, LinearAlgebra
import SparseArrays, LinearAlgebra
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Is an import necessary? Are we adding methods?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I looked at this the other way around:
import is safer than an unqualified using, with the only difference being that it avoids importing many unknown and unused methods into the namespace with possibility of name clashes etc?

In this case, the only effect of this change was to enforce SparseArrays.nnz which looks like what was intended, and caught a couple of places where an unqualified nnz was inadvertently used which then was potentially confusing (well confused me anyway!).

Many thanks for fantastic work making this package available btw!

Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

In this case, the only effect of this change was to enforce SparseArrays.nnz

Good point.

using ..SpkGraph: Graph
using ..SpkOrdering: Ordering
using ..SpkETree: ETree
Expand All @@ -8,7 +8,22 @@ using ..SpkSparseSolver: SparseSolver, findorder!,symbolicfactor!,inmatrix!,fact
import ..SpkSparseSolver: solve!
import ..SpkSparseBase: _inmatrix!

"""
Graph(m::SparseMatrixCSC, diagonal=false) -> g::Graph

Construct a graph from a problem object.

It does not check that the problem object contains a structurally
symmetric matrix, since sometimes only the lower or upper triangle of
a symmetric matrix may be stored. There are routines in the SpkGraph module to
make a given graph object structurally symmetric.

Input:
- `m` - the sparse matrix, used to create the graph
- `diagonal` - indicates that the diagonal elements are included. If
diagonal is not given, the adjacency structure does not include
the diagonal elements.
"""
function Graph(m::SparseMatrixCSC{FT,IT}, diagonal=false) where {FT,IT}
nv = size(m,1)
nrows = size(m,2)
Expand All @@ -18,7 +33,7 @@ function Graph(m::SparseMatrixCSC{FT,IT}, diagonal=false) where {FT,IT}


if (diagonal)
nedges = nnz(m)
nedges = SparseArrays.nnz(m)
else
dedges=0
for i in 1:ncols
Expand All @@ -29,7 +44,7 @@ function Graph(m::SparseMatrixCSC{FT,IT}, diagonal=false) where {FT,IT}
end
end
end
nedges = nnz(m) - dedges
nedges = SparseArrays.nnz(m) - dedges
end

#jf if diagonal == true, we possibly can just use colptr & rowval
Expand All @@ -56,8 +71,6 @@ function Graph(m::SparseMatrixCSC{FT,IT}, diagonal=false) where {FT,IT}
end




function _SparseBase(m::SparseMatrixCSC{FT,IT}) where {IT,FT}
maxblocksize = 30 # This can be set by the user

Expand Down Expand Up @@ -194,19 +207,19 @@ end
Solves linear system defined with sparse solver and provides the solution in rhs.
"""
function solve!(s::SparseSolver{IT,FT}, rhs) where {IT,FT}
findorder!(s) || ErrorException("Finding Order.")
symbolicfactor!(s) || ErrorException("Symbolic Factorization.")
inmatrix!(s) || ErrorException("Matrix input.")
factor!(s) || ErrorException("Numerical Factorization.")
triangularsolve!(s,rhs) || ErrorException("Triangular Solve.")
findorder!(s) || error("Finding Order.")
Copy link
Owner

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good catch with the exceptions!

symbolicfactor!(s) || error("Symbolic Factorization.")
inmatrix!(s) || error("Matrix input.")
factor!(s) || error("Numerical Factorization.")
triangularsolve!(s,rhs) || error("Triangular Solve.")
return true
end

#########################################################################
# SparspakLU

"""
sparspaklu(m; factorize=true)
sparspaklu(m::SparseMatrixCSC; factorize=true) -> lu::SparseSolver

If `factorize==true`, calculate LU factorization using Sparspak. Steps
are `findorder`, `symbolicfactor`, `factor`.
Expand All @@ -220,39 +233,52 @@ which has methods for `LinearAlgebra.ldiv!` and "backslash".
function sparspaklu(m::SparseMatrixCSC{FT,IT};factorize=true) where {FT,IT}
lu=SparseSolver(m)
if factorize
findorder!(lu) || ErrorException("Finding Order.")
symbolicfactor!(lu) || ErrorException("Symbolic Factorization.")
inmatrix!(lu) || ErrorException("Matrix input.")
factor!(lu) || ErrorException("Numerical Factorization.")
findorder!(lu) || error("Finding Order.")
symbolicfactor!(lu) || error("Symbolic Factorization.")
inmatrix!(lu) || error("Matrix input.")
factor!(lu) || error("Numerical Factorization.")
end
lu
end


"""
sparspaklu!(lu,m)
sparspaklu!(lu::SparseSolver, m::SparseMatrixCSC; allow_pattern_change=true) -> lu::SparseSolver

Calculate LU factorization of a sparse matrix `m`, reusing ordering and symbolic
factorization `lu`, if that was previously calculated.

Calculate LU factorization, reusing ordering and symbolic
factorization from lu, if that was previously calculated.
If `allow_pattern_change = true` (the default) the sparse matrix `m` may have a nonzero pattern
different to that of the matrix used to create the LU factorization `lu`, in which case the ordering
and symbolic factorization are updated.

Currently, it is assumed that, if size and number of nonzeros didn't
change, the sparsity patterns of `m` and `p` are the same, probably
leading to errors elsewhere if the patterns nevertheless differ.
If `allow_pattern_change = false` an error is thrown if the nonzero pattern of `m` is different to that
of the matrix used to create the LU factorization `lu`.

If `lu` has not been factorized (ie it has just been created with option `factorize = false`) then
`lu` is always updated from `m` and `allow_pattern_change` is ignored.

"""
function sparspaklu!(lu::SparseSolver{IT,FT}, m::SparseMatrixCSC{FT,IT}) where {FT,IT}
# jf: Do we need a better test here ? Not sure as that may be expensive.
if lu.slvr.n != size(m,1) || lu.slvr.n != size(m,2) || lu.slvr.nnz != nnz(m)
lu=SparseSolver(m)
end
function sparspaklu!(lu::SparseSolver{IT,FT}, m::SparseMatrixCSC{FT,IT}; allow_pattern_change=true) where {FT,IT}

pattern_changed = (SparseArrays.getcolptr(m) != SparseArrays.getcolptr(lu.p)) || (SparseArrays.getrowval(m) != SparseArrays.getrowval(lu.p))

if pattern_changed
if allow_pattern_change || !lu._symbolicdone
copy!(lu, SparseSolver(m))
else
error("'allow_pattern_change=false', but sparsity pattern of matrix 'm' does not match that used to create 'lu'")
end
end

lu.p=m
lu._orderdone || ( findorder!(lu) || ErrorException("Finding Order.") )
lu._symbolicdone || ( symbolicfactor!(lu) || ErrorException("Symbolic Factorization.") )
lu._orderdone || ( findorder!(lu) || error("Finding Order.") )
lu._symbolicdone || ( symbolicfactor!(lu) || error("Symbolic Factorization.") )
lu._inmatrixdone = false
lu._factordone = false
lu._trisolvedone = false
inmatrix!(lu) || ErrorException("Matrix input.")
factor!(lu) || ErrorException("Numerical Factorization.")
inmatrix!(lu) || error("Matrix input.")
factor!(lu) || error("Numerical Factorization.")
lu
end

Expand All @@ -276,7 +302,7 @@ Overwriting left division for SparseSolver.
The solution is returned in `v`, which is also the right hand side vector.
"""
function LinearAlgebra.ldiv!(lu::SparseSolver{IT,FT}, v) where {IT,FT}
solve!(lu, v) || ErrorException("Triangular Solve.")
solve!(lu, v) || error("Triangular Solve.")
lu._trisolvedone = false
v
end
Expand Down
17 changes: 12 additions & 5 deletions src/SparseMethod/SpkSparseSolver.jl
Original file line number Diff line number Diff line change
Expand Up @@ -31,6 +31,13 @@ mutable struct SparseSolver{IT, FT}
_condestdone::Bool
end

function Base.copy!(dst::SparseSolver, src::SparseSolver)
for f in fieldnames(SparseSolver)
setfield!(dst, f, getfield(src, f))
end
return dst
end

"""
SparseSolver(p::Problem)

Expand Down Expand Up @@ -71,11 +78,11 @@ Given a symmetric matrix `A`, the steps are:
The solution can be retrieved as `p.x`.
"""
function solve!(s::SparseSolver{IT}) where {IT}
findorder!(s) || ErrorException("Finding Order.")
symbolicfactor!(s) || ErrorException("Symbolic Factorization.")
inmatrix!(s) || ErrorException("Matrix input.")
factor!(s) || ErrorException("Numerical Factorization.")
triangularsolve!(s) || ErrorException("Triangular Solve.")
findorder!(s) || error("Finding Order.")
symbolicfactor!(s) || error("Symbolic Factorization.")
inmatrix!(s) || error("Matrix input.")
factor!(s) || error("Numerical Factorization.")
triangularsolve!(s) || error("Triangular Solve.")
return true
end

Expand Down
10 changes: 5 additions & 5 deletions src/SparseSpdMethod/SpkSparseSpdSolver.jl
Original file line number Diff line number Diff line change
Expand Up @@ -205,11 +205,11 @@ Given a symmetric matrix `A`, the steps are:
The solution can be retrieved as `p.x`.
"""
function solve!(s::SparseSpdSolver{IT}) where {IT}
findorder!(s) || ErrorException("Finding Order.")
symbolicfactor!(s) || ErrorException("Symbolic Factorization.")
inmatrix!(s) || ErrorException("Matrix input.")
factor!(s) || ErrorException("Numerical Factorization.")
triangularsolve!(s) || ErrorException("Triangular Solve.")
findorder!(s) || error("Finding Order.")
symbolicfactor!(s) || error("Symbolic Factorization.")
inmatrix!(s) || error("Matrix input.")
factor!(s) || error("Numerical Factorization.")
triangularsolve!(s) || error("Triangular Solve.")
return true
end

Expand Down
55 changes: 54 additions & 1 deletion test/test_cscinterface.jl
Original file line number Diff line number Diff line change
Expand Up @@ -137,7 +137,19 @@ function _test(T=Float64, n=20)

spm.nzval.-=0.1
rhs = spm*exsol
lu=sparspaklu!(lu,spm)
sparspaklu!(lu,spm)
sol=lu\rhs
@test sol≈exsol

# create a matrix with different sparsity pattern
spm2 = spm + sprand(T, n, n, 1/n)
rhs = spm2*exsol

# test fails attempting to reuse factorisation lu
@test_throws ErrorException sparspaklu!(lu,spm2; allow_pattern_change=false)

# test with default allow_pattern_change == true
sparspaklu!(lu,spm2)
sol=lu\rhs
@test sol≈exsol

Expand All @@ -147,4 +159,45 @@ end
_test(Float64)
_test(Float64x2)
_test(ForwardDiff.Dual{Float64,Float64,1})


function _test_asymmetric(T=Float64, n=20)
spm = sprand(T, n, n, 1/n)
spm = -spm + 40 * LinearAlgebra.I


exsol = ones(T,n)
rhs = spm*exsol
lu=sparspaklu(spm)
sol=lu\rhs
@test sol≈exsol

spm.nzval.-=0.1
rhs = spm*exsol
sparspaklu!(lu,spm)
sol=lu\rhs
@test sol≈exsol

# create a matrix with different sparsity pattern
spm2 = spm + sprand(T, n, n, 1/n)
rhs = spm2*exsol

# test fails attempting to reuse factorisation lu
@test_throws ErrorException sparspaklu!(lu,spm2; allow_pattern_change=false)

# test with default allow_pattern_change=true (will redo symbolic factorization)
sparspaklu!(lu,spm2)
sol=lu\rhs
@test sol≈exsol

# test with uninitialized lu
lu_nofact = sparspaklu(spzeros(1, 1); factorize=false)
sparspaklu!(lu_nofact,spm2; allow_pattern_change=false) # always allow update of unfactorized lu
sol=lu_nofact\rhs
@test sol≈exsol

end

_test_asymmetric(Float64)

end
Loading