Skip to content

Commit

Permalink
Don't throw in sparse LU for singular matrices. Fixes #13781.
Browse files Browse the repository at this point in the history
(cherry picked from commit 857a010)
ref #13784
  • Loading branch information
andreasnoack authored and tkelman committed Oct 31, 2015
1 parent f975d61 commit 25ab144
Show file tree
Hide file tree
Showing 2 changed files with 45 additions and 40 deletions.
82 changes: 42 additions & 40 deletions base/sparse/umfpack.jl
Original file line number Diff line number Diff line change
Expand Up @@ -12,45 +12,45 @@ import Base.SparseMatrix: increment, increment!, decrement, decrement!, nnz

include("umfpack_h.jl")
type MatrixIllConditionedException <: Exception
message :: AbstractString
message::AbstractString
end

function umferror(status::Integer)
if status==UMFPACK_OK
return
elseif status==UMFPACK_WARNING_singular_matrix
throw(MatrixIllConditionedException("singular matrix"))
elseif status==UMFPACK_WARNING_determinant_underflow
throw(MatrixIllConditionedException("the determinant is nonzero but underflowed"))
elseif status==UMFPACK_WARNING_determinant_overflow
throw(MatrixIllConditionedException("the determinant overflowed"))
elseif status==UMFPACK_ERROR_out_of_memory
throw(OutOfMemoryError())
elseif status==UMFPACK_ERROR_invalid_Numeric_object
throw(ArgumentError("invalid UMFPack numeric object"))
elseif status==UMFPACK_ERROR_invalid_Symbolic_object
throw(ArgumentError("invalid UMFPack symbolic object"))
elseif status==UMFPACK_ERROR_argument_missing
throw(ArgumentError("a required argument to UMFPack is missing"))
elseif status==UMFPACK_ERROR_n_nonpositive
throw(BoundsError("the number of rows or columns of the matrix must be greater than zero"))
elseif status==UMFPACK_ERROR_invalid_matrix
throw(ArgumentError("invalid matrix"))
elseif status==UMFPACK_ERROR_different_pattern
throw(ArgumentError("pattern of the matrix changed"))
elseif status==UMFPACK_ERROR_invalid_system
throw(ArgumentError("invalid sys argument provided to UMFPack solver"))
elseif status==UMFPACK_ERROR_invalid_permutation
throw(ArgumentError("invalid permutation"))
elseif status==UMFPACK_ERROR_file_IO
throw(ErrorException("error saving / loading UMFPack decomposition"))
elseif status==UMFPACK_ERROR_ordering_failed
throw(ErrorException("the ordering method failed"))
elseif status==UMFPACK_ERROR_internal_error
throw(ErrorException("an internal error has occurred, of unknown cause"))
else
throw(ErrorException("unknown UMFPack error code: $status"))
end
if status==UMFPACK_OK
return
elseif status==UMFPACK_WARNING_singular_matrix
throw(LinAlg.SingularException(0))
elseif status==UMFPACK_WARNING_determinant_underflow
throw(MatrixIllConditionedException("the determinant is nonzero but underflowed"))
elseif status==UMFPACK_WARNING_determinant_overflow
throw(MatrixIllConditionedException("the determinant overflowed"))
elseif status==UMFPACK_ERROR_out_of_memory
throw(OutOfMemoryError())
elseif status==UMFPACK_ERROR_invalid_Numeric_object
throw(ArgumentError("invalid UMFPack numeric object"))
elseif status==UMFPACK_ERROR_invalid_Symbolic_object
throw(ArgumentError("invalid UMFPack symbolic object"))
elseif status==UMFPACK_ERROR_argument_missing
throw(ArgumentError("a required argument to UMFPack is missing"))
elseif status==UMFPACK_ERROR_n_nonpositive
throw(BoundsError("the number of rows or columns of the matrix must be greater than zero"))
elseif status==UMFPACK_ERROR_invalid_matrix
throw(ArgumentError("invalid matrix"))
elseif status==UMFPACK_ERROR_different_pattern
throw(ArgumentError("pattern of the matrix changed"))
elseif status==UMFPACK_ERROR_invalid_system
throw(ArgumentError("invalid sys argument provided to UMFPack solver"))
elseif status==UMFPACK_ERROR_invalid_permutation
throw(ArgumentError("invalid permutation"))
elseif status==UMFPACK_ERROR_file_IO
throw(ErrorException("error saving / loading UMFPack decomposition"))
elseif status==UMFPACK_ERROR_ordering_failed
throw(ErrorException("the ordering method failed"))
elseif status==UMFPACK_ERROR_internal_error
throw(ErrorException("an internal error has occurred, of unknown cause"))
else
throw(ErrorException("unknown UMFPack error code: $status"))
end
end

macro isok(A)
Expand Down Expand Up @@ -170,8 +170,9 @@ for itype in UmfpackIndexTypes
Ptr{Float64}, Ptr{Float64}),
U.colptr, U.rowval, U.nzval, U.symbolic, tmp,
umf_ctrl, umf_info)
status > 0 && throw(MatrixIllConditionedException(""))
umferror(status)
if status != UMFPACK_WARNING_singular_matrix
umferror(status)
end
U.numeric = tmp[1]
return U
end
Expand All @@ -184,8 +185,9 @@ for itype in UmfpackIndexTypes
Ptr{Float64}, Ptr{Float64}),
U.colptr, U.rowval, real(U.nzval), imag(U.nzval), U.symbolic, tmp,
umf_ctrl, umf_info)
status > 0 && throw(MatrixIllConditionedException(""))
umferror(status)
if status != UMFPACK_WARNING_singular_matrix
umferror(status)
end
U.numeric = tmp[1]
return U
end
Expand Down
3 changes: 3 additions & 0 deletions test/sparsedir/umfpack.jl
Original file line number Diff line number Diff line change
Expand Up @@ -59,3 +59,6 @@ end
#4523 - complex sparse \
x = speye(2) + im * speye(2)
@test_approx_eq (x*(lufact(x) \ ones(2))) ones(2)

@test det(sparse([1,3,3,1], [1,1,3,3], [1,1,1,1])) == 0

0 comments on commit 25ab144

Please sign in to comment.