diff --git a/src/factorization.jl b/src/factorization.jl index 0ef7b754e..efcaaaf00 100644 --- a/src/factorization.jl +++ b/src/factorization.jl @@ -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 """ @@ -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 diff --git a/test/default_algs.jl b/test/default_algs.jl index d2aafba9d..c604129c9 100644 --- a/test/default_algs.jl +++ b/test/default_algs.jl @@ -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))