Skip to content

Commit 988e3f7

Browse files
committed
Use Pivoted QR for Underdetermined Systems
1 parent 86beacf commit 988e3f7

File tree

4 files changed

+46
-5
lines changed

4 files changed

+46
-5
lines changed

src/LinearSolve.jl

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -62,6 +62,9 @@ needs_concrete_A(alg::AbstractKrylovSubspaceMethod) = false
6262
needs_concrete_A(alg::AbstractSolveFunction) = false
6363

6464
# Util
65+
is_underdetermined(x) = false
66+
is_underdetermined(A::AbstractMatrix) = size(A, 1) < size(A, 2)
67+
is_underdetermined(A::AbstractSciMLOperator) = size(A, 1) < size(A, 2)
6568

6669
_isidentity_struct(A) = false
6770
_isidentity_struct::Number) = isone(λ)
@@ -96,6 +99,7 @@ EnumX.@enumx DefaultAlgorithmChoice begin
9699
NormalCholeskyFactorization
97100
AppleAccelerateLUFactorization
98101
MKLLUFactorization
102+
QRFactorizationPivoted
99103
end
100104

101105
struct DefaultLinearSolver <: SciMLLinearSolveAlgorithm

src/default.jl

Lines changed: 23 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
needs_concrete_A(alg::DefaultLinearSolver) = true
22
mutable struct DefaultLinearSolverInit{T1, T2, T3, T4, T5, T6, T7, T8, T9, T10, T11, T12,
3-
T13, T14, T15, T16, T17, T18}
3+
T13, T14, T15, T16, T17, T18, T19}
44
LUFactorization::T1
55
QRFactorization::T2
66
DiagonalFactorization::T3
@@ -19,6 +19,7 @@ mutable struct DefaultLinearSolverInit{T1, T2, T3, T4, T5, T6, T7, T8, T9, T10,
1919
NormalCholeskyFactorization::T16
2020
AppleAccelerateLUFactorization::T17
2121
MKLLUFactorization::T18
22+
QRFactorizationPivoted::T19
2223
end
2324

2425
# Legacy fallback
@@ -168,8 +169,8 @@ function defaultalg(A, b, assump::OperatorAssumptions)
168169
(A === nothing ? eltype(b) <: Union{Float32, Float64} :
169170
eltype(A) <: Union{Float32, Float64})
170171
DefaultAlgorithmChoice.RFLUFactorization
171-
#elseif A === nothing || A isa Matrix
172-
# alg = FastLUFactorization()
172+
#elseif A === nothing || A isa Matrix
173+
# alg = FastLUFactorization()
173174
elseif usemkl && (A === nothing ? eltype(b) <: Union{Float32, Float64} :
174175
eltype(A) <: Union{Float32, Float64})
175176
DefaultAlgorithmChoice.MKLLUFactorization
@@ -199,9 +200,19 @@ function defaultalg(A, b, assump::OperatorAssumptions)
199200
elseif assump.condition === OperatorCondition.WellConditioned
200201
DefaultAlgorithmChoice.NormalCholeskyFactorization
201202
elseif assump.condition === OperatorCondition.IllConditioned
202-
DefaultAlgorithmChoice.QRFactorization
203+
if is_underdetermined(A)
204+
# Underdetermined
205+
DefaultAlgorithmChoice.QRFactorizationPivoted
206+
else
207+
DefaultAlgorithmChoice.QRFactorization
208+
end
203209
elseif assump.condition === OperatorCondition.VeryIllConditioned
204-
DefaultAlgorithmChoice.QRFactorization
210+
if is_underdetermined(A)
211+
# Underdetermined
212+
DefaultAlgorithmChoice.QRFactorizationPivoted
213+
else
214+
DefaultAlgorithmChoice.QRFactorization
215+
end
205216
elseif assump.condition === OperatorCondition.SuperIllConditioned
206217
DefaultAlgorithmChoice.SVDFactorization
207218
else
@@ -247,6 +258,12 @@ function algchoice_to_alg(alg::Symbol)
247258
NormalCholeskyFactorization()
248259
elseif alg === :AppleAccelerateLUFactorization
249260
AppleAccelerateLUFactorization()
261+
elseif alg === :QRFactorizationPivoted
262+
@static if VERSION v"1.7beta"
263+
QRFactorization(ColumnNorm())
264+
else
265+
QRFactorization(Val(true))
266+
end
250267
else
251268
error("Algorithm choice symbol $alg not allowed in the default")
252269
end
@@ -310,6 +327,7 @@ function defaultalg_symbol(::Type{T}) where {T}
310327
Symbol(split(string(SciMLBase.parameterless_type(T)), ".")[end])
311328
end
312329
defaultalg_symbol(::Type{<:GenericFactorization{typeof(ldlt!)}}) = :LDLtFactorization
330+
defaultalg_symbol(::Type{<:QRFactorization{ColumnNorm}}) = :QRFactorizationPivoted
313331

314332
"""
315333
if alg.alg === DefaultAlgorithmChoice.LUFactorization

src/factorization.jl

Lines changed: 10 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -158,6 +158,16 @@ function QRFactorization(inplace = true)
158158
QRFactorization(pivot, 16, inplace)
159159
end
160160

161+
@static if VERSION v"1.7beta"
162+
function QRFactorization(pivot::LinearAlgebra.PivotingStrategy, inplace::Bool = true)
163+
QRFactorization(pivot, 16, inplace)
164+
end
165+
else
166+
function QRFactorization(pivot::Val, inplace::Bool = true)
167+
QRFactorization(pivot, 16, inplace)
168+
end
169+
end
170+
161171
function do_factorization(alg::QRFactorization, A, b, u)
162172
A = convert(AbstractMatrix, A)
163173
if alg.inplace && !(A isa SparseMatrixCSC) && !(A isa GPUArraysCore.AbstractGPUArray)

test/nonsquare.jl

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -56,3 +56,12 @@ solve(LinearProblem(A, b), (LinearSolve.NormalCholeskyFactorization())).u;
5656
solve(LinearProblem(A, b),
5757
assumptions = (OperatorAssumptions(false;
5858
condition = OperatorCondition.WellConditioned))).u;
59+
60+
# Underdetermined
61+
m, n = 2, 3
62+
63+
A = rand(m, n)
64+
b = rand(m)
65+
prob = LinearProblem(A, b)
66+
res = A \ b
67+
@test solve(prob).u res

0 commit comments

Comments
 (0)