Skip to content

Commit 725325c

Browse files
Fix and better test UMFPACK default
Fixes #326
1 parent b5ff8b1 commit 725325c

File tree

2 files changed

+20
-3
lines changed

2 files changed

+20
-3
lines changed

src/factorization.jl

Lines changed: 4 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -707,16 +707,17 @@ function SciMLBase.solve!(cache::LinearCache, alg::UMFPACKFactorization; kwargs.
707707
A = cache.A
708708
A = convert(AbstractMatrix, A)
709709
if cache.isfresh
710+
cacheval = @get_cacheval(cache, :UMFPACKFactorization)
710711
if alg.reuse_symbolic
711712
# Caches the symbolic factorization: https://github.com/JuliaLang/julia/pull/33738
712713
if alg.check_pattern && !(SuiteSparse.decrement(SparseArrays.getcolptr(A)) ==
713-
cache.cacheval.colptr &&
714+
cacheval.colptr &&
714715
SuiteSparse.decrement(SparseArrays.getrowval(A)) ==
715-
@get_cacheval(cache, :UMFPACKFactorization).rowval)
716+
cacheval.rowval)
716717
fact = lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
717718
nonzeros(A)))
718719
else
719-
fact = lu!(@get_cacheval(cache, :UMFPACKFactorization),
720+
fact = lu!(cacheval,
720721
SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
721722
nonzeros(A)))
722723
end

test/default_algs.jl

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,10 +1,19 @@
11
using LinearSolve, LinearAlgebra, SparseArrays, Test, JET
22
@test LinearSolve.defaultalg(nothing, zeros(3)).alg ===
33
LinearSolve.DefaultAlgorithmChoice.GenericLUFactorization
4+
prob = LinearProblem(rand(3,3), rand(3))
5+
solve(prob)
6+
47
@test LinearSolve.defaultalg(nothing, zeros(50)).alg ===
58
LinearSolve.DefaultAlgorithmChoice.RFLUFactorization
9+
prob = LinearProblem(rand(50,50), rand(50))
10+
solve(prob)
11+
612
@test LinearSolve.defaultalg(nothing, zeros(600)).alg ===
713
LinearSolve.DefaultAlgorithmChoice.GenericLUFactorization
14+
prob = LinearProblem(rand(600,600), rand(600))
15+
solve(prob)
16+
817
@test LinearSolve.defaultalg(LinearAlgebra.Diagonal(zeros(5)), zeros(5)).alg ===
918
LinearSolve.DefaultAlgorithmChoice.DiagonalFactorization
1019

@@ -14,8 +23,15 @@ using LinearSolve, LinearAlgebra, SparseArrays, Test, JET
1423

1524
@test LinearSolve.defaultalg(sprand(1000, 1000, 0.01), zeros(1000)).alg ===
1625
LinearSolve.DefaultAlgorithmChoice.KLUFactorization
26+
prob = LinearProblem(sprand(1000, 1000, 0.01), zeros(1000))
27+
solve(prob)
28+
1729
@test LinearSolve.defaultalg(sprand(11000, 11000, 0.001), zeros(11000)).alg ===
1830
LinearSolve.DefaultAlgorithmChoice.UMFPACKFactorization
31+
prob = LinearProblem(sprand(11000, 11000, 0.001), zeros(11000))
32+
solve(prob)
33+
34+
1935

2036
@static if VERSION >= v"v1.7-"
2137
# Test inference

0 commit comments

Comments
 (0)