Skip to content

Commit cbc6178

Browse files
Merge pull request #75 from JuliaLinearAlgebra/rms-lulinv-with-pivoting
Include pivoting in LUL⁻¹
2 parents d30ad01 + 7275b86 commit cbc6178

File tree

3 files changed

+131
-74
lines changed

3 files changed

+131
-74
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "MatrixFactorizations"
22
uuid = "a3b82374-2e81-5b9e-98ce-41277c0e4c87"
3-
version = "3.1.1"
3+
version = "3.1.2"
44

55
[deps]
66
ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a"

src/lulinv.jl

Lines changed: 77 additions & 44 deletions
Original file line numberDiff line numberDiff line change
@@ -10,6 +10,8 @@ The individual components of the factorization `F::LULinv` can be accessed via [
1010
|:----------|:--------------------------------------------|
1111
| `F.L` | `L` (unit lower triangular) part of `LUL⁻¹` |
1212
| `F.U` | `U` (upper triangular) part of `LUL⁻¹` |
13+
| `F.p` | (right) permutation `Vector` |
14+
| `F.P` | (right) permutation `Matrix` |
1315
1416
Iterating the factorization produces the components `F.L` and `F.U`.
1517
@@ -31,12 +33,12 @@ U factor:
3133
-0.772002 3.0
3234
0.0 7.772
3335
34-
julia> F.L * F.U / F.L ≈ A
36+
julia> F.L * F.U / F.L ≈ A[F.p, F.p]
3537
true
3638
37-
julia> l, u = lulinv(A); # destructuring via iteration
39+
julia> l, u, p = lulinv(A); # destructuring via iteration
3840
39-
julia> l == F.L && u == F.U
41+
julia> l == F.L && u == F.U && p == F.p
4042
true
4143
4244
julia> A = [-150 334 778; -89 195 464; 5 -10 -27]
@@ -58,25 +60,27 @@ U factor:
5860
0 -2 75
5961
0 0 3
6062
61-
julia> F.L * F.U / F.L == A
63+
julia> F.L * F.U / F.L == A[F.p, F.p]
6264
true
6365
```
6466
"""
65-
struct LULinv{T, S <: AbstractMatrix{T}} <: Factorization{T}
67+
struct LULinv{T, S <: AbstractMatrix{T}, IPIV <: AbstractVector{<:Integer}} <: Factorization{T}
6668
factors::S
67-
function LULinv{T, S}(factors) where {T, S <: AbstractMatrix{T}}
69+
ipiv::IPIV
70+
function LULinv{T, S, IPIV}(factors, ipiv) where {T, S <: AbstractMatrix{T}, IPIV <: AbstractVector{<:Integer}}
6871
require_one_based_indexing(factors)
69-
new{T, S}(factors)
72+
new{T, S, IPIV}(factors, ipiv)
7073
end
7174
end
7275

7376

74-
LULinv(factors::AbstractMatrix{T}) where T = LULinv{T, typeof(factors)}(factors)
75-
LULinv{T}(factors::AbstractMatrix) where T = LULinv(convert(AbstractMatrix{T}, factors))
76-
LULinv{T}(F::LULinv) where T = LULinv{T}(F.factors)
77+
LULinv(factors::AbstractMatrix{T}, ipiv::AbstractVector{<:Integer}) where T = LULinv{T, typeof(factors), typeof(ipiv)}(factors, ipiv)
78+
LULinv{T}(factors::AbstractMatrix, ipiv::AbstractVector{<:Integer}) where T = LULinv(convert(AbstractMatrix{T}, factors), ipiv)
79+
LULinv{T}(F::LULinv) where T = LULinv{T}(F.factors, F.ipiv)
7780

7881
iterate(F::LULinv) = (F.L, Val(:U))
79-
iterate(F::LULinv, ::Val{:U}) = (F.U, Val(:done))
82+
iterate(F::LULinv, ::Val{:U}) = (F.U, Val(:p))
83+
iterate(F::LULinv, ::Val{:p}) = (F.p, Val(:done))
8084
iterate(F::LULinv, ::Val{:done}) = nothing
8185

8286

@@ -116,17 +120,22 @@ function getU(F::LULinv)
116120
end
117121

118122
function getproperty(F::LULinv{T, <: AbstractMatrix}, d::Symbol) where T
123+
n = size(F, 1)
119124
if d === :L
120125
return getL(F)
121126
elseif d === :U
122127
return getU(F)
128+
elseif d === :p
129+
return getfield(F, :ipiv)
130+
elseif d === :P
131+
return Matrix{T}(I, n, n)[:, F.p]
123132
else
124133
getfield(F, d)
125134
end
126135
end
127136

128137
propertynames(F::LULinv, private::Bool=false) =
129-
(:L, :U, (private ? fieldnames(typeof(F)) : ())...)
138+
(:L, :U, :p, :P, (private ? fieldnames(typeof(F)) : ())...)
130139

131140
function show(io::IO, mime::MIME{Symbol("text/plain")}, F::LULinv)
132141
summary(io, F); println(io)
@@ -137,42 +146,60 @@ function show(io::IO, mime::MIME{Symbol("text/plain")}, F::LULinv)
137146
end
138147

139148

140-
lulinv(A::AbstractMatrix; kwds...) = lulinv(A, eigvals(A); kwds...)
141-
function lulinv(A::AbstractMatrix{T}, λ::AbstractVector{T}; kwds...) where T
149+
lulinv(A::AbstractMatrix, pivot::Union{Val{false}, Val{true}}=Val(true); kwds...) = lulinv(A, eigvals(A), pivot; kwds...)
150+
function lulinv(A::AbstractMatrix{T}, λ::AbstractVector{T}, pivot::Union{Val{false}, Val{true}}=Val(true); kwds...) where T
142151
S = lulinvtype(T)
143-
lulinv!(copy_oftype(A, S), copy_oftype(λ, S); kwds...)
152+
lulinv!(copy_oftype(A, S), copy_oftype(λ, S), pivot; kwds...)
144153
end
145-
function lulinv(A::AbstractMatrix{T1}, λ::AbstractVector{T2}; kwds...) where {T1, T2}
154+
function lulinv(A::AbstractMatrix{T1}, λ::AbstractVector{T2}, pivot::Union{Val{false}, Val{true}}=Val(true); kwds...) where {T1, T2}
146155
T = promote_type(T1, T2)
147156
S = lulinvtype(T)
148-
lulinv!(copy_oftype(A, S), copy_oftype(λ, S); kwds...)
157+
lulinv!(copy_oftype(A, S), copy_oftype(λ, S), pivot; kwds...)
149158
end
150159

151-
function lulinv!(A::Matrix{T}, λ::Vector{T}; rtol::Real = size(A, 1)*eps(real(float(oneunit(T))))) where T
160+
function lulinv!(A::Matrix{T}, λ::Vector{T}, ::Val{Pivot} = Val(true); rtol::Real = size(A, 1)*eps(real(float(oneunit(T))))) where {T, Pivot}
152161
n = checksquare(A)
153162
n == length(λ) || throw(ArgumentError("Eigenvalue count does not match matrix dimensions."))
154163
v = zeros(T, n)
155-
for i in 1:n-1
156-
# We must find an eigenvector with nonzero "first" entry. A failed UL factorization of A-λI reveals this vector provided L₁₁ == 0.
157-
for j in 1:length(λ)
158-
F = ul!(view(A, i:n, i:n) - λ[j]*I; check=false)
159-
nrm = norm(F.L)
160-
if norm(F.L[1]) rtol*nrm # v[i] is a free parameter; we set it to 1.
161-
fill!(v, zero(T))
162-
v[i] = one(T)
163-
# Next, we must scan the remaining free parameters, set them to 0, so that we find a nonsingular lower-triangular linear system for the nontrivial remaining part of the eigenvector.
164-
idx = Int[]
165-
for k in 2:n+1-i
166-
if norm(F.L[k, k]) > rtol*nrm
167-
push!(idx, k)
168-
end
164+
ipiv = collect(1:n)
165+
for i in 1:n
166+
F = ul!(view(A, i:n, i:n) - λ[i]*I; check = false)
167+
nrm = norm(F.L)
168+
if Pivot
169+
ip = n+1-i
170+
while (ip > 1) && (norm(F.L[ip, ip]) > rtol*nrm)
171+
ip -= 1
172+
end
173+
fill!(v, zero(T))
174+
v[i] = one(T)
175+
idx = ip+1:n+1-i
176+
if !isempty(idx)
177+
v[idx.+(i-1)] .= -F.L[idx, ip]
178+
ldiv!(LowerTriangular(view(F.L, idx, idx)), view(v, idx.+(i-1)))
179+
end
180+
ip = ip+i-1
181+
if ip i
182+
tmp = ipiv[i]
183+
ipiv[i] = ipiv[ip]
184+
ipiv[ip] = tmp
185+
for j in 1:n
186+
tmp = A[i, j]
187+
A[i, j] = A[ip, j]
188+
A[ip, j] = tmp
169189
end
170-
if !isempty(idx)
171-
v[idx.+(i-1)] .= -F.L[idx, 1]
172-
ldiv!(LowerTriangular(view(F.L, idx, idx)), view(v, idx.+(i-1)))
190+
for j in 1:n
191+
tmp = A[j, i]
192+
A[j, i] = A[j, ip]
193+
A[j, ip] = tmp
173194
end
174-
deleteat!(λ, j)
175-
break
195+
end
196+
else
197+
fill!(v, zero(T))
198+
v[i] = one(T)
199+
idx = 2:n+1-i
200+
if !isempty(idx)
201+
v[idx.+(i-1)] .= -F.L[idx, 1]
202+
ldiv!(LowerTriangular(view(F.L, idx, idx)), view(v, idx.+(i-1)))
176203
end
177204
end
178205
for k in 1:n
@@ -189,17 +216,23 @@ function lulinv!(A::Matrix{T}, λ::Vector{T}; rtol::Real = size(A, 1)*eps(real(f
189216
A[k, i] = v[k]
190217
end
191218
end
192-
return LULinv(A)
219+
return LULinv(A, ipiv)
193220
end
194221

195-
function ldiv!(F::LULinv, B::AbstractVecOrMat)
196-
L = UnitLowerTriangular(F.factors)
197-
return lmul!(L, ldiv!(UpperTriangular(F.factors), ldiv!(L, B)))
222+
function ldiv!(A::LULinv, B::AbstractVecOrMat)
223+
B .= B[invperm(A.p), :] # Todo: fix me
224+
L = UnitLowerTriangular(A.factors)
225+
lmul!(L, ldiv!(UpperTriangular(A.factors), ldiv!(L, B)))
226+
B .= B[A.p, :] # Todo: fix me
227+
return B
198228
end
199229

200-
function rdiv!(B::AbstractVecOrMat, F::LULinv)
201-
L = UnitLowerTriangular(F.factors)
202-
return rdiv!(rdiv!(rmul!(B, L), UpperTriangular(F.factors)), L)
230+
function rdiv!(A::AbstractVecOrMat, B::LULinv)
231+
A .= A[:, B.p] # Todo: fix me
232+
L = UnitLowerTriangular(B.factors)
233+
rdiv!(rdiv!(rmul!(A, L), UpperTriangular(B.factors)), L)
234+
A .= A[:, invperm(B.p)] # Todo: fix me
235+
return A
203236
end
204237

205238
det(F::LULinv) = det(UpperTriangular(F.factors))

test/test_lulinv.jl

Lines changed: 53 additions & 29 deletions
Original file line numberDiff line numberDiff line change
@@ -6,18 +6,18 @@ using LinearAlgebra, MatrixFactorizations, Random, Test
66
λ = [17//1; -2; 3]
77
A = det(V)*V*Diagonal(λ)/V
88
@test A == Matrix{Int}(A)
9-
L, U = lulinv(A, λ)
10-
@test A L*U/L
11-
@test A*L L*U
9+
L, U, p = lulinv(A, λ)
10+
@test A[p, p] L*U/L
11+
@test A[p, p]*L L*U
1212
end
1313
@testset "A::Matrix{Rational{Int}}" begin
1414
Random.seed!(0)
1515
V = rand(-3:3//1, 5, 5)
1616
λ = [-2; -1; 0; 1; 2//1]
1717
A = V*Diagonal(λ)/V
18-
L, U = lulinv(A, λ)
19-
@test A L*U/L
20-
@test A*L L*U
18+
L, U, p = lulinv(A, λ)
19+
@test A[p, p] L*U/L
20+
@test A[p, p]*L L*U
2121
end
2222
@testset "Full Jordan blocks" begin
2323
Random.seed!(0)
@@ -28,9 +28,10 @@ using LinearAlgebra, MatrixFactorizations, Random, Test
2828
J[3, 4] = 1
2929
J[4, 5] = 1
3030
A = V*J/V
31-
L, U = lulinv(A, λ)
32-
@test A L*U/L
33-
@test A*L L*U
31+
L, U, p = lulinv(A, λ)
32+
@test A[p, p] L*U/L
33+
@test A[p, p]*L L*U
34+
3435
end
3536
@testset "Alg ≠ geo" begin
3637
Random.seed!(0)
@@ -40,17 +41,17 @@ using LinearAlgebra, MatrixFactorizations, Random, Test
4041
J[1, 2] = 1
4142
J[3, 4] = 1
4243
A = V*J/V
43-
L, U = lulinv(A, λ)
44-
@test A L*U/L
45-
@test A*L L*U
44+
L, U, p = lulinv(A, λ)
45+
@test A[p, p] L*U/L
46+
@test A[p, p]*L L*U
4647
λ = [1; -2; -2; 1; 1//1]
47-
L, U = lulinv(A, λ)
48-
@test A L*U/L
49-
@test A*L L*U
48+
L, U, p = lulinv(A, λ)
49+
@test A[p, p] L*U/L
50+
@test A[p, p]*L L*U
5051
λ = [1; -2; 1; -2; 1//1]
51-
L, U = lulinv(A, λ)
52-
@test A L*U/L
53-
@test A*L L*U
52+
L, U, p = lulinv(A, λ)
53+
@test A[p, p] L*U/L
54+
@test A[p, p]*L L*U
5455
end
5556
@testset "Index scan is correct" begin
5657
Random.seed!(0)
@@ -60,15 +61,15 @@ using LinearAlgebra, MatrixFactorizations, Random, Test
6061
J = [J1 zeros(Rational{Int}, size(J1)); zeros(Rational{Int}, size(J2)) J2]
6162
λ = diag(J)
6263
A = V*J/V
63-
L, U = lulinv(A, λ)
64-
@test A L*U/L
65-
@test A*L L*U
64+
L, U, p = lulinv(A, λ)
65+
@test A[p, p] L*U/L
66+
@test A[p, p]*L L*U
6667
end
6768
@testset "A::Matrix{Float64}" begin
6869
A = [4.0 3; 6 3]
69-
L, U = lulinv(A)
70-
@test A L*U/L
71-
@test A*L L*U
70+
L, U, p = lulinv(A)
71+
@test A[p, p] L*U/L
72+
@test A[p, p]*L L*U
7273
end
7374
@testset "size, det, logdet, logabsdet" begin
7475
A = [4.0 3; 6 3]
@@ -90,7 +91,7 @@ using LinearAlgebra, MatrixFactorizations, Random, Test
9091
seekstart(bf)
9192
@test String(take!(bf)) ==
9293
"""
93-
LULinv{Float64, Matrix{Float64}}
94+
LULinv{Float64, Matrix{Float64}, Vector{Int64}}
9495
L factor:
9596
2×2 Matrix{Float64}:
9697
1.0 0.0
@@ -102,15 +103,38 @@ U factor:
102103
end
103104
@testset "propertynames" begin
104105
names = sort!(collect(string.(Base.propertynames(lulinv([2 1; 1 2])))))
105-
@test names == ["L", "U"]
106+
@test names == ["L", "P", "U", "p"]
106107
allnames = sort!(collect(string.(Base.propertynames(lulinv([2 1; 1 2]), true))))
107-
@test allnames == ["L", "U", "factors"]
108+
@test allnames == ["L", "P", "U", "factors", "ipiv", "p"]
108109
end
109110
@testset "A::Matrix{Int64}, λ::Vector{Rational{Int64}}" begin
110111
A = [-150 334 778; -89 195 464; 5 -10 -27]
111112
F = lulinv(A, [17, -2, 3//1])
112-
@test A * F.L == F.L * F.U
113+
@test A[F.p, F.p] * F.L == F.L * F.U
113114
G = LULinv{Float64}(F)
114-
@test A * G.L G.L * G.U
115+
@test A[G.p, G.p] * G.L G.L * G.U
116+
end
117+
@testset "A is lower triangular" begin
118+
Random.seed!(0)
119+
for A in (tril!(rand(-5:5, 8, 8)), tril!(randn(8, 8)), tril!(randn(Float32, 8, 8)))
120+
F = lulinv(A)
121+
@test A[F.p, F.p] * F.L F.L * F.U
122+
L, U, P = F.L, F.U, F.P
123+
@test P'A*P*L L*U
124+
@test P'A*P L*U/L
125+
@test A P*L*U/(P*L)
126+
@test A P*L*U/L*P'
127+
end
128+
end
129+
@testset "Pivoting is necessary" begin
130+
A = [1 0; 0 2//1]
131+
λ = [2; 1//1]
132+
@test_throws SingularException lulinv(A, λ, Val(false))
133+
L, U, p = lulinv(A, λ)
134+
@test p == [2, 1]
135+
@test A[p, p] == L*U/L
136+
F = lulinv(A, λ)
137+
@test F\A I
138+
@test A/F I
115139
end
116140
end

0 commit comments

Comments
 (0)