@@ -6,7 +6,7 @@ export UmfpackLU
6
6
7
7
import Base: (\ ), getproperty, show, size
8
8
using LinearAlgebra
9
- import LinearAlgebra: Factorization, det, lu, ldiv!
9
+ import LinearAlgebra: Factorization, det, lu, lu!, ldiv!
10
10
11
11
using SparseArrays
12
12
using SparseArrays: getcolptr
@@ -176,6 +176,65 @@ lu(A::Union{SparseMatrixCSC{T},SparseMatrixCSC{Complex{T}}};
176
176
" dense LU." )))
177
177
lu (A:: SparseMatrixCSC ; check:: Bool = true ) = lu (float (A); check = check)
178
178
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)
179
238
180
239
size (F:: UmfpackLU ) = (F. m, F. n)
181
240
function size (F:: UmfpackLU , dim:: Integer )
@@ -262,35 +321,33 @@ for itype in UmfpackIndexTypes
262
321
return U
263
322
end
264
323
function umfpack_numeric! (U:: UmfpackLU{Float64,$itype} )
265
- if U. numeric != C_NULL return U end
266
324
if U. symbolic == C_NULL umfpack_symbolic! (U) end
267
325
tmp = Vector {Ptr{Cvoid}} (undef, 1 )
268
326
status = ccall (($ num_r, :libumfpack ), $ itype,
269
327
(Ptr{$ itype}, Ptr{$ itype}, Ptr{Float64}, Ptr{Cvoid}, Ptr{Cvoid},
270
328
Ptr{Float64}, Ptr{Float64}),
271
329
U. colptr, U. rowval, U. nzval, U. symbolic, tmp,
272
330
umf_ctrl, umf_info)
331
+ U. status = status
273
332
if status != UMFPACK_WARNING_singular_matrix
274
333
umferror (status)
275
334
end
276
335
U. numeric = tmp[1 ]
277
- U. status = status
278
336
return U
279
337
end
280
338
function umfpack_numeric! (U:: UmfpackLU{ComplexF64,$itype} )
281
- if U. numeric != C_NULL return U end
282
339
if U. symbolic == C_NULL umfpack_symbolic! (U) end
283
340
tmp = Vector {Ptr{Cvoid}} (undef, 1 )
284
341
status = ccall (($ num_c, :libumfpack ), $ itype,
285
342
(Ptr{$ itype}, Ptr{$ itype}, Ptr{Float64}, Ptr{Float64}, Ptr{Cvoid}, Ptr{Cvoid},
286
343
Ptr{Float64}, Ptr{Float64}),
287
344
U. colptr, U. rowval, real (U. nzval), imag (U. nzval), U. symbolic, tmp,
288
345
umf_ctrl, umf_info)
346
+ U. status = status
289
347
if status != UMFPACK_WARNING_singular_matrix
290
348
umferror (status)
291
349
end
292
350
U. numeric = tmp[1 ]
293
- U. status = status
294
351
return U
295
352
end
296
353
function solve! (x:: StridedVector{Float64} , lu:: UmfpackLU{Float64,$itype} , b:: StridedVector{Float64} , typ:: Integer )
0 commit comments