Skip to content

Commit

Permalink
Merge pull request #474 from oscardssmith/remove-klu-and-umfpack-from…
Browse files Browse the repository at this point in the history
…-default

Fix `UMFPACK` to use `check=false`
  • Loading branch information
ChrisRackauckas authored Feb 24, 2024
2 parents e23f5c7 + 86e7191 commit fc5c14a
Show file tree
Hide file tree
Showing 2 changed files with 19 additions and 6 deletions.
17 changes: 11 additions & 6 deletions src/factorization.jl
Original file line number Diff line number Diff line change
Expand Up @@ -780,21 +780,26 @@ function SciMLBase.solve!(cache::LinearCache, alg::UMFPACKFactorization; kwargs.
SparseArrays.decrement(SparseArrays.getrowval(A)) ==
cacheval.rowval)
fact = lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)))
nonzeros(A)), check=false)
else
fact = lu!(cacheval,
SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)))
nonzeros(A)), check=false)
end
else
fact = lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A), nonzeros(A)))
fact = lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A), nonzeros(A)), check=false)
end
cache.cacheval = fact
cache.isfresh = false
end

y = ldiv!(cache.u, @get_cacheval(cache, :UMFPACKFactorization), cache.b)
SciMLBase.build_linear_solution(alg, y, nothing, cache)
F = @get_cacheval(cache, :UMFPACKFactorization)
if F.status == SparseArrays.UMFPACK.UMFPACK_OK
y = ldiv!(cache.u, F, cache.b)
SciMLBase.build_linear_solution(alg, y, nothing, cache)
else
SciMLBase.build_linear_solution(alg, cache.u, nothing, cache; retcode=ReturnCode.Infeasible)
end
end

"""
Expand Down Expand Up @@ -840,10 +845,10 @@ function init_cacheval(alg::KLUFactorization, A::AbstractSparseArray, b, u, Pl,
nonzeros(A)))
end

# TODO: guard this against errors
function SciMLBase.solve!(cache::LinearCache, alg::KLUFactorization; kwargs...)
A = cache.A
A = convert(AbstractMatrix, A)

if cache.isfresh
cacheval = @get_cacheval(cache, :KLUFactorization)
if cacheval !== nothing && alg.reuse_symbolic
Expand Down
8 changes: 8 additions & 0 deletions test/default_algs.jl
Original file line number Diff line number Diff line change
Expand Up @@ -35,6 +35,14 @@ solve(prob)
LinearSolve.OperatorAssumptions(false)).alg ===
LinearSolve.DefaultAlgorithmChoice.QRFactorization


A = spzeros(2, 2)
# test that solving a singular problem doesn't error
prob = LinearProblem(A, ones(2))
@test solve(prob, UMFPACKFactorization()).retcode == ReturnCode.Infeasible
@test_broken solve(prob, KLUFactorization()).retcode == ReturnCode.Infeasible


@test LinearSolve.defaultalg(sprand(10^4, 10^4, 1e-5) + I, zeros(1000)).alg ===
LinearSolve.DefaultAlgorithmChoice.KLUFactorization
prob = LinearProblem(sprand(1000, 1000, 0.5), zeros(1000))
Expand Down

0 comments on commit fc5c14a

Please sign in to comment.