Skip to content

Commit

Permalink
Add check=false with newest KLU (#478)
Browse files Browse the repository at this point in the history
* Update to most recent KLU

* Unbreak test

* fix accidental dev leak
  • Loading branch information
rayegun authored Feb 29, 2024
1 parent f97cf52 commit 6d51feb
Show file tree
Hide file tree
Showing 3 changed files with 18 additions and 19 deletions.
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ HYPRE = "1.4.0"
InteractiveUtils = "1.10"
IterativeSolvers = "0.9.3"
JET = "0.8.28"
KLU = "0.5"
KLU = "0.6"
KernelAbstractions = "0.9.16"
Krylov = "0.9"
KrylovKit = "0.6"
Expand Down
33 changes: 16 additions & 17 deletions src/factorization.jl
Original file line number Diff line number Diff line change
Expand Up @@ -855,23 +855,17 @@ function SciMLBase.solve!(cache::LinearCache, alg::KLUFactorization; kwargs...)
A = convert(AbstractMatrix, A)
if cache.isfresh
cacheval = @get_cacheval(cache, :KLUFactorization)
if cacheval !== nothing && alg.reuse_symbolic
if alg.reuse_symbolic
if alg.check_pattern && !(SparseArrays.decrement(SparseArrays.getcolptr(A)) ==
cacheval.colptr &&
SparseArrays.decrement(SparseArrays.getrowval(A)) == cacheval.rowval)
fact = KLU.klu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)))
else
# If we have a cacheval already, run umfpack_symbolic to ensure the symbolic factorization exists
# This won't recompute if it does.
KLU.klu_analyze!(cacheval)
copyto!(cacheval.nzval, nonzeros(A))
if cacheval._numeric === C_NULL # We MUST have a numeric factorization for reuse, unlike UMFPACK.
KLU.klu_factor!(cacheval)
end
fact = KLU.klu!(cacheval,
SparseArrays.decrement(SparseArrays.getrowval(A)) ==
cacheval.rowval)
fact = KLU.klu(
SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)))
nonzeros(A)),
check = false)
else
fact = KLU.klu!(cacheval, nonzeros(A), check = false)
end
else
# New fact each time since the sparsity pattern can change
Expand All @@ -882,9 +876,14 @@ function SciMLBase.solve!(cache::LinearCache, alg::KLUFactorization; kwargs...)
cache.cacheval = fact
cache.isfresh = false
end

y = ldiv!(cache.u, @get_cacheval(cache, :KLUFactorization), cache.b)
SciMLBase.build_linear_solution(alg, y, nothing, cache)
F = @get_cacheval(cache, :KLUFactorization)
if F.common.status == KLU.KLU_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

## CHOLMODFactorization
Expand Down
2 changes: 1 addition & 1 deletion test/default_algs.jl
Original file line number Diff line number Diff line change
Expand Up @@ -39,7 +39,7 @@ 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 solve(prob, KLUFactorization()).retcode == ReturnCode.Infeasible

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

0 comments on commit 6d51feb

Please sign in to comment.