Skip to content

Commit c31e460

Browse files
Merge pull request #253 from j-fu/update-sparspak
Update sparspak
2 parents 7ad09fc + bae7cc0 commit c31e460

File tree

5 files changed

+81
-18
lines changed

5 files changed

+81
-18
lines changed

Project.toml

Lines changed: 4 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -39,15 +39,17 @@ Reexport = "1"
3939
SciMLBase = "1.68"
4040
Setfield = "0.7, 0.8, 1"
4141
SnoopPrecompile = "1"
42-
Sparspak = "0.3"
42+
Sparspak = "0.3.4"
4343
UnPack = "1"
4444
julia = "1.6"
4545

4646
[extras]
47+
MultiFloats = "bdf0d083-296b-4888-a5b6-7498122e68a5"
48+
ForwardDiff = "f6369f11-7733-5829-9624-2563aa707210"
4749
Pkg = "44cfe95a-1eb2-52ea-b672-e2afdf69b78f"
4850
Random = "9a3f8284-a2c9-5f02-9a11-845980a1fd5c"
4951
SafeTestsets = "1bc83da4-3b8d-516f-aca4-4fe02f6d838f"
5052
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
5153

5254
[targets]
53-
test = ["Test", "Pkg", "Random", "SafeTestsets"]
55+
test = ["Test", "Pkg", "Random", "SafeTestsets", "MultiFloats", "ForwardDiff"]

docs/src/solvers/solvers.md

Lines changed: 20 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -21,6 +21,12 @@ For sparse LU-factorizations, `KLUFactorization` if there is less structure
2121
to the sparsity pattern and `UMFPACKFactorization` if there is more structure.
2222
Pardiso.jl's methods are also known to be very efficient sparse linear solvers.
2323

24+
While these sparse factorizations are based on implementations in other languages,
25+
and therefore constrained to standard number types (`Float64`, `Float32` and
26+
their complex counterparts), `SparspakFactorization` is able to handle general
27+
number types, e.g. defined by `ForwardDiff.jl`, `MultiFloats.jl`,
28+
or `IntervalArithmetics.jl`.
29+
2430
As sparse matrices get larger, iterative solvers tend to get more efficient than
2531
factorization methods if a lower tolerance of the solution is required.
2632

@@ -147,6 +153,20 @@ Base.@kwdef struct PardisoJL <: SciMLLinearSolveAlgorithm
147153
end
148154
```
149155

156+
### Sparspak.jl
157+
This is the translation of the well-known sparse matrix software Sparspak
158+
(Waterloo Sparse Matrix Package), solving
159+
large sparse systems of linear algebraic equations. Sparspak is composed of the
160+
subroutines from the book "Computer Solution of Large Sparse Positive Definite
161+
Systems" by Alan George and Joseph Liu. Originally written in Fortran 77, later
162+
rewritten in Fortran 90. Here is the software translated into Julia.
163+
The Julia rewrite is released under the MIT license with an express permission
164+
from the authors of the Fortran package. The package uses mutiple
165+
dispatch to route around standard BLAS routines in the case e.g. of arbitrary-precision
166+
floating point numbers or ForwardDiff.Dual.
167+
This e.g. allows for Automatic Differentiation (AD) of a sparse-matrix solve.
168+
169+
150170
### CUDA.jl
151171

152172
Note that `CuArrays` are supported by `GenericFactorization` in the “normal” way.

src/LinearSolve.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -78,6 +78,7 @@ end
7878
sol = solve(prob)
7979
sol = solve(prob, KLUFactorization())
8080
sol = solve(prob, UMFPACKFactorization())
81+
sol = solve(prob, SparspakFactorization())
8182
end
8283
end
8384

src/factorization.jl

Lines changed: 12 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -503,30 +503,28 @@ end
503503

504504
## SparspakFactorization is here since it's MIT licensed, not GPL
505505

506-
struct SparspakFactorization <: AbstractFactorization end
506+
Base.@kwdef struct SparspakFactorization <: AbstractFactorization
507+
reuse_symbolic::Bool = true
508+
end
507509

508510
function init_cacheval(::SparspakFactorization, A, b, u, Pl, Pr, maxiters::Int, abstol,
509511
reltol,
510512
verbose::Bool, assumptions::OperatorAssumptions)
511513
A = convert(AbstractMatrix, A)
512-
p = Sparspak.Problem.Problem(size(A)...)
513-
Sparspak.Problem.insparse!(p, A)
514-
Sparspak.Problem.infullrhs!(p, b)
515-
s = Sparspak.SparseSolver.SparseSolver(p)
516-
return s
514+
return sparspaklu(A, factorize=false)
517515
end
518516

519517
function SciMLBase.solve(cache::LinearCache, alg::SparspakFactorization; kwargs...)
520518
A = cache.A
521519
A = convert(AbstractMatrix, A)
522520
if cache.isfresh
523-
p = Sparspak.Problem.Problem(size(A)...)
524-
Sparspak.Problem.insparse!(p, A)
525-
Sparspak.Problem.infullrhs!(p, cache.b)
526-
s = Sparspak.SparseSolver.SparseSolver(p)
527-
cache = set_cacheval(cache, s)
521+
if cache.cacheval !== nothing && alg.reuse_symbolic
522+
fact = sparspaklu!(cache.cacheval, A)
523+
else
524+
fact = sparspaklu(A)
525+
end
526+
cache = set_cacheval(cache, fact)
528527
end
529-
Sparspak.SparseSolver.solve!(cache.cacheval)
530-
copyto!(cache.u, cache.cacheval.p.x)
531-
SciMLBase.build_linear_solution(alg, cache.u, nothing, cache)
528+
y = ldiv!(cache.u, cache.cacheval, cache.b)
529+
SciMLBase.build_linear_solution(alg, y, nothing, cache)
532530
end

test/basictests.jl

Lines changed: 44 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,9 @@
1-
using LinearSolve, LinearAlgebra, SparseArrays
1+
using LinearSolve, LinearAlgebra, SparseArrays, MultiFloats, ForwardDiff
22
using Test
33
import Random
44

5+
const Dual64=ForwardDiff.Dual{Nothing,Float64,1}
6+
57
n = 8
68
A = Matrix(I, n, n)
79
b = ones(n)
@@ -124,7 +126,7 @@ end
124126
@test X * solve(cache) b1
125127
end
126128

127-
@testset "Sparspak Factorization" begin
129+
@testset "Sparspak Factorization (Float64)" begin
128130
A1 = sparse(A / 1)
129131
b1 = rand(n)
130132
x1 = zero(b)
@@ -137,6 +139,46 @@ end
137139
test_interface(SparspakFactorization(), prob1, prob2)
138140
end
139141

142+
@testset "Sparspak Factorization (Float64x1)" begin
143+
A1 = sparse(A / 1) .|> Float64x1
144+
b1 = rand(n) .|> Float64x1
145+
x1 = zero(b) .|> Float64x1
146+
A2 = sparse(A / 2) .|> Float64x1
147+
b2 = rand(n) .|> Float64x1
148+
x2 = zero(b) .|> Float64x1
149+
150+
prob1 = LinearProblem(A1, b1; u0 = x1)
151+
prob2 = LinearProblem(A2, b2; u0 = x2)
152+
test_interface(SparspakFactorization(), prob1, prob2)
153+
end
154+
155+
@testset "Sparspak Factorization (Float64x2)" begin
156+
A1 = sparse(A / 1) .|> Float64x2
157+
b1 = rand(n) .|> Float64x2
158+
x1 = zero(b) .|> Float64x2
159+
A2 = sparse(A / 2) .|> Float64x2
160+
b2 = rand(n) .|> Float64x2
161+
x2 = zero(b) .|> Float64x2
162+
163+
prob1 = LinearProblem(A1, b1; u0 = x1)
164+
prob2 = LinearProblem(A2, b2; u0 = x2)
165+
test_interface(SparspakFactorization(), prob1, prob2)
166+
end
167+
168+
@testset "Sparspak Factorization (Dual64)" begin
169+
A1 = sparse(A / 1) .|> Dual64
170+
b1 = rand(n) .|> Dual64
171+
x1 = zero(b) .|> Dual64
172+
A2 = sparse(A / 2) .|> Dual64
173+
b2 = rand(n) .|> Dual64
174+
x2 = zero(b) .|> Dual64
175+
176+
prob1 = LinearProblem(A1, b1; u0 = x1)
177+
prob2 = LinearProblem(A2, b2; u0 = x2)
178+
test_interface(SparspakFactorization(), prob1, prob2)
179+
end
180+
181+
140182
@testset "FastLAPACK Factorizations" begin
141183
A1 = A / 1
142184
b1 = rand(n)

0 commit comments

Comments
 (0)