Skip to content
Open
24 changes: 19 additions & 5 deletions src/factorization.jl
Original file line number Diff line number Diff line change
Expand Up @@ -723,8 +723,14 @@ end

################################## Factorizations which require solve! overloads

# Default control array of 20 elements from Umfpack

#!!!ATTENTION
#These array could change and whenever there will be a change in the SparseArrays package this array must be changed
const default_control=[1.0, 0.2, 0.2, 0.1, 32.0, 0.0, 0.7, 0.0, 1.0, 0.3, 1.0, 1.0, 0.9, 0.0, 10.0, 0.001, 1.0, 0.5, 0.0, 1.0]

"""
`UMFPACKFactorization(;reuse_symbolic=true, check_pattern=true)`
`UMFPACKFactorization(;reuse_symbolic=true, check_pattern=true, control=Vector{Float64}(20))`

A fast sparse multithreaded LU-factorization which specializes on sparsity
patterns with “more structure”.
Expand All @@ -735,10 +741,17 @@ patterns with “more structure”.
symbolic factorization. I.e., if `set_A` is used, it is expected that the new
`A` has the same sparsity pattern as the previous `A`. If this algorithm is to
be used in a context where that assumption does not hold, set `reuse_symbolic=false`.

There is also the possibility to set your own control array for the factorization that will be
passed at the SparseArrays package. To do so u have to pass a Vector{Float64} with 20 elements.
If no vecetor is passed the default control array is used.
More details on the UMFPACK doc
"""
Base.@kwdef struct UMFPACKFactorization <: AbstractFactorization

Base.@kwdef struct UMFPACKFactorization{T <: Union{Nothing, Vector{Float64}}} <: AbstractFactorization
reuse_symbolic::Bool = true
check_pattern::Bool = true # Check factorization re-use
control::T=nothing
end

const PREALLOCATED_UMFPACK = SparseArrays.UMFPACK.UmfpackLU(SparseMatrixCSC(0, 0, [1],
Expand Down Expand Up @@ -771,6 +784,7 @@ end
function SciMLBase.solve!(cache::LinearCache, alg::UMFPACKFactorization; kwargs...)
A = cache.A
A = convert(AbstractMatrix, A)
isnothing(alg.control) ? control=default_control : control=alg.control
if cache.isfresh
cacheval = @get_cacheval(cache, :UMFPACKFactorization)
if alg.reuse_symbolic
Expand All @@ -782,15 +796,15 @@ function SciMLBase.solve!(cache::LinearCache, alg::UMFPACKFactorization; kwargs.
fact = lu(
SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)),
check = false)
check = false, control=control)
else
fact = lu!(cacheval,
SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A),
nonzeros(A)), check = false)
nonzeros(A)), check = false, control=control)
end
else
fact = lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A), nonzeros(A)),
check = false)
check = false, control=control)
end
cache.cacheval = fact
cache.isfresh = false
Expand Down
12 changes: 12 additions & 0 deletions test/controlvectortest.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,12 @@
using SparseArrays, LinearSolve

A = sparse(rand(3,3));
b=rand(3);
prob = LinearProblem(A, b);

#check without control Vector
u=solve(prob,UMFPACKFactorization()).u

#check plugging in a control vector
controlv=SparseArrays.UMFPACK.get_umfpack_control(Float64,Int64)
u=solve(prob,UMFPACKFactorization(control=controlv)).u