Skip to content

Commit 1eb4dee

Browse files
goggleViralBShah
authored andcommitted
Introduce lu! for UmfpackLU to make use of its symbolic factorization (#33738)
* Introduce `lu!` for `UmfpackLU` to make use of its symbolic factorization * Remove spaces at keyword argument Co-authored-by: Viral B. Shah <viral@juliacomputing.com>
1 parent cc18589 commit 1eb4dee

File tree

3 files changed

+108
-5
lines changed

3 files changed

+108
-5
lines changed

NEWS.md

Lines changed: 18 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -43,6 +43,24 @@ Standard library changes
4343

4444

4545
#### SparseArrays
46+
* `lu!` accepts `UmfpackLU` as an argument to make use of its symbolic factorization.
47+
48+
#### Dates
49+
50+
#### Statistics
51+
52+
53+
#### Sockets
54+
55+
56+
Deprecated or removed
57+
---------------------
58+
59+
External dependencies
60+
---------------------
61+
62+
Tooling Improvements
63+
---------------------
4664

4765

4866
<!--- generated by NEWS-update.jl: -->

stdlib/SuiteSparse/src/umfpack.jl

Lines changed: 62 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -6,7 +6,7 @@ export UmfpackLU
66

77
import Base: (\), getproperty, show, size
88
using LinearAlgebra
9-
import LinearAlgebra: Factorization, det, lu, ldiv!
9+
import LinearAlgebra: Factorization, det, lu, lu!, ldiv!
1010

1111
using SparseArrays
1212
using SparseArrays: getcolptr
@@ -176,6 +176,65 @@ lu(A::Union{SparseMatrixCSC{T},SparseMatrixCSC{Complex{T}}};
176176
"dense LU.")))
177177
lu(A::SparseMatrixCSC; check::Bool = true) = lu(float(A); check = check)
178178

179+
"""
180+
lu!(F::UmfpackLU, A::SparseMatrixCSC; check=true) -> F::UmfpackLU
181+
182+
Compute the LU factorization of a sparse matrix `A`, reusing the symbolic
183+
factorization of an already existing LU factorization stored in `F`. The
184+
sparse matrix `A` must have an identical nonzero pattern as the matrix used
185+
to create the LU factorization `F`, otherwise an error is thrown.
186+
187+
When `check = true`, an error is thrown if the decomposition fails.
188+
When `check = false`, responsibility for checking the decomposition's
189+
validity (via [`issuccess`](@ref)) lies with the user.
190+
191+
!!! note
192+
`lu!(F::UmfpackLU, A::SparseMatrixCSC)` uses the UMFPACK library that is part of
193+
SuiteSparse. As this library only supports sparse matrices with [`Float64`](@ref) or
194+
`ComplexF64` elements, `lu!` converts `A` into a copy that is of type
195+
`SparseMatrixCSC{Float64}` or `SparseMatrixCSC{ComplexF64}` as appropriate.
196+
197+
!!! compat "Julia 1.4"
198+
`lu!` for `UmfpackLU` requires at least Julia 1.4.
199+
200+
# Examples
201+
```jldoctest
202+
julia> A = sparse(Float64[1.0 2.0; 0.0 3.0]);
203+
204+
julia> F = lu(A);
205+
206+
julia> B = sparse(Float64[1.0 1.0; 0.0 1.0]);
207+
208+
julia> lu!(F, B);
209+
210+
julia> F \\ ones(2)
211+
2-element Array{Float64,1}:
212+
0.0
213+
1.0
214+
```
215+
"""
216+
function lu!(F::UmfpackLU, S::SparseMatrixCSC{<:UMFVTypes,<:UMFITypes}; check::Bool=true)
217+
zerobased = getcolptr(S)[1] == 0
218+
F.m = size(S, 1)
219+
F.n = size(S, 2)
220+
F.colptr = zerobased ? copy(getcolptr(S)) : decrement(getcolptr(S))
221+
F.rowval = zerobased ? copy(rowvals(S)) : decrement(rowvals(S))
222+
F.nzval = copy(nonzeros(S))
223+
224+
umfpack_numeric!(F)
225+
check && (issuccess(F) || throw(LinearAlgebra.SingularException(0)))
226+
return F
227+
end
228+
lu!(F::UmfpackLU, A::SparseMatrixCSC{<:Union{Float16,Float32},Ti};
229+
check::Bool = true) where {Ti<:UMFITypes} =
230+
lu!(F, convert(SparseMatrixCSC{Float64,Ti}, A); check = check)
231+
lu!(F::UmfpackLU, A::SparseMatrixCSC{<:Union{ComplexF16,ComplexF32},Ti};
232+
check::Bool = true) where {Ti<:UMFITypes} =
233+
lu!(F, convert(SparseMatrixCSC{ComplexF64,Ti}, A); check = check)
234+
lu!(F::UmfpackLU, A::Union{SparseMatrixCSC{T},SparseMatrixCSC{Complex{T}}};
235+
check::Bool = true) where {T<:AbstractFloat} =
236+
throw(ArgumentError(string("matrix type ", typeof(A), "not supported.")))
237+
lu!(F::UmfpackLU, A::SparseMatrixCSC; check::Bool = true) = lu!(F, float(A); check = check)
179238

180239
size(F::UmfpackLU) = (F.m, F.n)
181240
function size(F::UmfpackLU, dim::Integer)
@@ -262,35 +321,33 @@ for itype in UmfpackIndexTypes
262321
return U
263322
end
264323
function umfpack_numeric!(U::UmfpackLU{Float64,$itype})
265-
if U.numeric != C_NULL return U end
266324
if U.symbolic == C_NULL umfpack_symbolic!(U) end
267325
tmp = Vector{Ptr{Cvoid}}(undef, 1)
268326
status = ccall(($num_r, :libumfpack), $itype,
269327
(Ptr{$itype}, Ptr{$itype}, Ptr{Float64}, Ptr{Cvoid}, Ptr{Cvoid},
270328
Ptr{Float64}, Ptr{Float64}),
271329
U.colptr, U.rowval, U.nzval, U.symbolic, tmp,
272330
umf_ctrl, umf_info)
331+
U.status = status
273332
if status != UMFPACK_WARNING_singular_matrix
274333
umferror(status)
275334
end
276335
U.numeric = tmp[1]
277-
U.status = status
278336
return U
279337
end
280338
function umfpack_numeric!(U::UmfpackLU{ComplexF64,$itype})
281-
if U.numeric != C_NULL return U end
282339
if U.symbolic == C_NULL umfpack_symbolic!(U) end
283340
tmp = Vector{Ptr{Cvoid}}(undef, 1)
284341
status = ccall(($num_c, :libumfpack), $itype,
285342
(Ptr{$itype}, Ptr{$itype}, Ptr{Float64}, Ptr{Float64}, Ptr{Cvoid}, Ptr{Cvoid},
286343
Ptr{Float64}, Ptr{Float64}),
287344
U.colptr, U.rowval, real(U.nzval), imag(U.nzval), U.symbolic, tmp,
288345
umf_ctrl, umf_info)
346+
U.status = status
289347
if status != UMFPACK_WARNING_singular_matrix
290348
umferror(status)
291349
end
292350
U.numeric = tmp[1]
293-
U.status = status
294351
return U
295352
end
296353
function solve!(x::StridedVector{Float64}, lu::UmfpackLU{Float64,$itype}, b::StridedVector{Float64}, typ::Integer)

stdlib/SuiteSparse/test/umfpack.jl

Lines changed: 28 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -189,6 +189,34 @@ using LinearAlgebra: Adjoint, Transpose, SingularException
189189
end
190190
end
191191

192+
@testset "Reuse symbolic LU factorization" begin
193+
A1 = sparse(increment!([0,4,1,1,2,2,0,1,2,3,4,4]),
194+
increment!([0,4,0,2,1,2,1,4,3,2,1,2]),
195+
[2.,1.,3.,4.,-1.,-3.,3.,9.,2.,1.,4.,2.], 5, 5)
196+
for Tv in (Float64, ComplexF64, Float16, Float32, ComplexF16, ComplexF32)
197+
for Ti in Base.uniontypes(SuiteSparse.UMFPACK.UMFITypes)
198+
A = convert(SparseMatrixCSC{Tv,Ti}, A0)
199+
B = convert(SparseMatrixCSC{Tv,Ti}, A1)
200+
b = Tv[8., 45., -3., 3., 19.]
201+
F = lu(A)
202+
lu!(F, B)
203+
@test F\b B\b Matrix(B)\b
204+
205+
# singular matrix
206+
C = copy(B)
207+
C[4, 3] = Tv(0)
208+
F = lu(A)
209+
@test_throws SingularException lu!(F, C)
210+
211+
# change of nonzero pattern
212+
D = copy(B)
213+
D[5, 1] = Tv(1.0)
214+
F = lu(A)
215+
@test_throws ArgumentError lu!(F, D)
216+
end
217+
end
218+
end
219+
192220
end
193221

194222
@testset "REPL printing of UmfpackLU" begin

0 commit comments

Comments
 (0)