/
bunchkaufman.jl
394 lines (351 loc) · 12.6 KB
/
bunchkaufman.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
# This file is a part of Julia. License is MIT: https://julialang.org/license
## Create an extractor that extracts the modified original matrix, e.g.
## LD for BunchKaufman, UL for CholeskyDense, LU for LUDense and
## define size methods for Factorization types using it.
"""
BunchKaufman <: Factorization
Matrix factorization type of the Bunch-Kaufman factorization of a symmetric or
Hermitian matrix `A` as `P'UDU'P` or `P'LDL'P`, depending on whether the upper
(the default) or the lower triangle is stored in `A`. If `A` is complex symmetric
then `U'` and `L'` denote the unconjugated transposes, i.e. `transpose(U)` and
`transpose(L)`, respectively. This is the return type of [`bunchkaufman`](@ref),
the corresponding matrix factorization function.
If `S::BunchKaufman` is the factorization object, the components can be obtained
via `S.D`, `S.U` or `S.L` as appropriate given `S.uplo`, and `S.p`.
Iterating the decomposition produces the components `S.D`, `S.U` or `S.L`
as appropriate given `S.uplo`, and `S.p`.
# Examples
```jldoctest
julia> A = [1 2; 2 3]
2×2 Matrix{Int64}:
1 2
2 3
julia> S = bunchkaufman(A) # A gets wrapped internally by Symmetric(A)
BunchKaufman{Float64, Matrix{Float64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
-0.333333 0.0
0.0 3.0
U factor:
2×2 UnitUpperTriangular{Float64, Matrix{Float64}}:
1.0 0.666667
⋅ 1.0
permutation:
2-element Vector{Int64}:
1
2
julia> d, u, p = S; # destructuring via iteration
julia> d == S.D && u == S.U && p == S.p
true
julia> S = bunchkaufman(Symmetric(A, :L))
BunchKaufman{Float64, Matrix{Float64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
3.0 0.0
0.0 -0.333333
L factor:
2×2 UnitLowerTriangular{Float64, Matrix{Float64}}:
1.0 ⋅
0.666667 1.0
permutation:
2-element Vector{Int64}:
2
1
```
"""
struct BunchKaufman{T,S<:AbstractMatrix} <: Factorization{T}
LD::S
ipiv::Vector{BlasInt}
uplo::Char
symmetric::Bool
rook::Bool
info::BlasInt
function BunchKaufman{T,S}(LD, ipiv, uplo, symmetric, rook, info) where {T,S<:AbstractMatrix}
require_one_based_indexing(LD)
new(LD, ipiv, uplo, symmetric, rook, info)
end
end
BunchKaufman(A::AbstractMatrix{T}, ipiv::Vector{BlasInt}, uplo::AbstractChar, symmetric::Bool,
rook::Bool, info::BlasInt) where {T} =
BunchKaufman{T,typeof(A)}(A, ipiv, uplo, symmetric, rook, info)
# iteration for destructuring into components
Base.iterate(S::BunchKaufman) = (S.D, Val(:UL))
Base.iterate(S::BunchKaufman, ::Val{:UL}) = (S.uplo == 'L' ? S.L : S.U, Val(:p))
Base.iterate(S::BunchKaufman, ::Val{:p}) = (S.p, Val(:done))
Base.iterate(S::BunchKaufman, ::Val{:done}) = nothing
"""
bunchkaufman!(A, rook::Bool=false; check = true) -> BunchKaufman
`bunchkaufman!` is the same as [`bunchkaufman`](@ref), but saves space by overwriting the
input `A`, instead of creating a copy.
"""
function bunchkaufman!(A::RealHermSymComplexSym{T,S} where {T<:BlasReal,S<:StridedMatrix},
rook::Bool = false; check::Bool = true)
LD, ipiv, info = rook ? LAPACK.sytrf_rook!(A.uplo, A.data) : LAPACK.sytrf!(A.uplo, A.data)
check && checknonsingular(info)
BunchKaufman(LD, ipiv, A.uplo, true, rook, info)
end
function bunchkaufman!(A::Hermitian{T,S} where {T<:BlasComplex,S<:StridedMatrix{T}},
rook::Bool = false; check::Bool = true)
LD, ipiv, info = rook ? LAPACK.hetrf_rook!(A.uplo, A.data) : LAPACK.hetrf!(A.uplo, A.data)
check && checknonsingular(info)
BunchKaufman(LD, ipiv, A.uplo, false, rook, info)
end
function bunchkaufman!(A::StridedMatrix{<:BlasFloat}, rook::Bool = false; check::Bool = true)
if ishermitian(A)
return bunchkaufman!(Hermitian(A), rook; check = check)
elseif issymmetric(A)
return bunchkaufman!(Symmetric(A), rook; check = check)
else
throw(ArgumentError("Bunch-Kaufman decomposition is only valid for symmetric or Hermitian matrices"))
end
end
"""
bunchkaufman(A, rook::Bool=false; check = true) -> S::BunchKaufman
Compute the Bunch-Kaufman [^Bunch1977] factorization of a symmetric or
Hermitian matrix `A` as `P'*U*D*U'*P` or `P'*L*D*L'*P`, depending on
which triangle is stored in `A`, and return a [`BunchKaufman`](@ref) object.
Note that if `A` is complex symmetric then `U'` and `L'` denote
the unconjugated transposes, i.e. `transpose(U)` and `transpose(L)`.
Iterating the decomposition produces the components `S.D`, `S.U` or `S.L`
as appropriate given `S.uplo`, and `S.p`.
If `rook` is `true`, rook pivoting is used. If `rook` is false,
rook pivoting is not used.
When `check = true`, an error is thrown if the decomposition fails.
When `check = false`, responsibility for checking the decomposition's
validity (via [`issuccess`](@ref)) lies with the user.
The following functions are available for `BunchKaufman` objects:
[`size`](@ref), `\\`, [`inv`](@ref), [`issymmetric`](@ref),
[`ishermitian`](@ref), [`getindex`](@ref).
[^Bunch1977]: J R Bunch and L Kaufman, Some stable methods for calculating inertia and solving symmetric linear systems, Mathematics of Computation 31:137 (1977), 163-179. [url](http://www.ams.org/journals/mcom/1977-31-137/S0025-5718-1977-0428694-0/).
# Examples
```jldoctest
julia> A = [1 2; 2 3]
2×2 Matrix{Int64}:
1 2
2 3
julia> S = bunchkaufman(A) # A gets wrapped internally by Symmetric(A)
BunchKaufman{Float64, Matrix{Float64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
-0.333333 0.0
0.0 3.0
U factor:
2×2 UnitUpperTriangular{Float64, Matrix{Float64}}:
1.0 0.666667
⋅ 1.0
permutation:
2-element Vector{Int64}:
1
2
julia> d, u, p = S; # destructuring via iteration
julia> d == S.D && u == S.U && p == S.p
true
julia> S.U*S.D*S.U' - S.P*A*S.P'
2×2 Matrix{Float64}:
0.0 0.0
0.0 0.0
julia> S = bunchkaufman(Symmetric(A, :L))
BunchKaufman{Float64, Matrix{Float64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
3.0 0.0
0.0 -0.333333
L factor:
2×2 UnitLowerTriangular{Float64, Matrix{Float64}}:
1.0 ⋅
0.666667 1.0
permutation:
2-element Vector{Int64}:
2
1
julia> S.L*S.D*S.L' - A[S.p, S.p]
2×2 Matrix{Float64}:
0.0 0.0
0.0 0.0
```
"""
bunchkaufman(A::AbstractMatrix{T}, rook::Bool=false; check::Bool = true) where {T} =
bunchkaufman!(copy_oftype(A, typeof(sqrt(oneunit(T)))), rook; check = check)
BunchKaufman{T}(B::BunchKaufman) where {T} =
BunchKaufman(convert(Matrix{T}, B.LD), B.ipiv, B.uplo, B.symmetric, B.rook, B.info)
Factorization{T}(B::BunchKaufman) where {T} = BunchKaufman{T}(B)
size(B::BunchKaufman) = size(getfield(B, :LD))
size(B::BunchKaufman, d::Integer) = size(getfield(B, :LD), d)
issymmetric(B::BunchKaufman) = B.symmetric
ishermitian(B::BunchKaufman{T}) where T = T<:Real || !B.symmetric
function _ipiv2perm_bk(v::AbstractVector{T}, maxi::Integer, uplo::AbstractChar, rook::Bool) where T
require_one_based_indexing(v)
p = T[1:maxi;]
uploL = uplo == 'L'
i = uploL ? 1 : maxi
# if uplo == 'U' we construct the permutation backwards
@inbounds while 1 <= i <= length(v)
vi = v[i]
if vi > 0 # the 1x1 blocks
p[i], p[vi] = p[vi], p[i]
i += uploL ? 1 : -1
else # the 2x2 blocks
if rook
p[i], p[-vi] = p[-vi], p[i]
end
if uploL
vp = rook ? -v[i+1] : -vi
p[i + 1], p[vp] = p[vp], p[i + 1]
i += 2
else # 'U'
vp = rook ? -v[i-1] : -vi
p[i - 1], p[vp] = p[vp], p[i - 1]
i -= 2
end
end
end
return p
end
function getproperty(B::BunchKaufman{T}, d::Symbol) where {T<:BlasFloat}
n = size(B, 1)
if d === :p
return _ipiv2perm_bk(getfield(B, :ipiv), n, getfield(B, :uplo), B.rook)
elseif d === :P
return Matrix{T}(I, n, n)[:,invperm(B.p)]
elseif d === :L || d === :U || d === :D
if getfield(B, :rook)
LUD, od = LAPACK.syconvf_rook!(getfield(B, :uplo), 'C', copy(getfield(B, :LD)), getfield(B, :ipiv))
else
LUD, od = LAPACK.syconv!(getfield(B, :uplo), copy(getfield(B, :LD)), getfield(B, :ipiv))
end
if d === :D
if getfield(B, :uplo) == 'L'
odl = od[1:n - 1]
return Tridiagonal(odl, diag(LUD), getfield(B, :symmetric) ? odl : conj.(odl))
else # 'U'
odu = od[2:n]
return Tridiagonal(getfield(B, :symmetric) ? odu : conj.(odu), diag(LUD), odu)
end
elseif d === :L
if getfield(B, :uplo) == 'L'
return UnitLowerTriangular(LUD)
else
throw(ArgumentError("factorization is U*D*transpose(U) but you requested L"))
end
else # :U
if B.uplo == 'U'
return UnitUpperTriangular(LUD)
else
throw(ArgumentError("factorization is L*D*transpose(L) but you requested U"))
end
end
else
getfield(B, d)
end
end
Base.propertynames(B::BunchKaufman, private::Bool=false) =
(:p, :P, :L, :U, :D, (private ? fieldnames(typeof(B)) : ())...)
issuccess(B::BunchKaufman) = B.info == 0
function adjoint(B::BunchKaufman)
if ishermitian(B)
return B
else
throw(ArgumentError("adjoint not implemented for complex symmetric matrices"))
end
end
function Base.show(io::IO, mime::MIME{Symbol("text/plain")}, B::BunchKaufman)
if issuccess(B)
summary(io, B); println(io)
println(io, "D factor:")
show(io, mime, B.D)
println(io, "\n$(B.uplo) factor:")
show(io, mime, B.uplo == 'L' ? B.L : B.U)
println(io, "\npermutation:")
show(io, mime, B.p)
else
print(io, "Failed factorization of type $(typeof(B))")
end
end
function inv(B::BunchKaufman{<:BlasReal})
if B.rook
copytri!(LAPACK.sytri_rook!(B.uplo, copy(B.LD), B.ipiv), B.uplo, true)
else
copytri!(LAPACK.sytri!(B.uplo, copy(B.LD), B.ipiv), B.uplo, true)
end
end
function inv(B::BunchKaufman{<:BlasComplex})
if issymmetric(B)
if B.rook
copytri!(LAPACK.sytri_rook!(B.uplo, copy(B.LD), B.ipiv), B.uplo)
else
copytri!(LAPACK.sytri!(B.uplo, copy(B.LD), B.ipiv), B.uplo)
end
else
if B.rook
copytri!(LAPACK.hetri_rook!(B.uplo, copy(B.LD), B.ipiv), B.uplo, true)
else
copytri!(LAPACK.hetri!(B.uplo, copy(B.LD), B.ipiv), B.uplo, true)
end
end
end
function ldiv!(B::BunchKaufman{T}, R::StridedVecOrMat{T}) where T<:BlasReal
if B.rook
LAPACK.sytrs_rook!(B.uplo, B.LD, B.ipiv, R)
else
LAPACK.sytrs!(B.uplo, B.LD, B.ipiv, R)
end
end
function ldiv!(B::BunchKaufman{T}, R::StridedVecOrMat{T}) where T<:BlasComplex
if B.rook
if issymmetric(B)
LAPACK.sytrs_rook!(B.uplo, B.LD, B.ipiv, R)
else
LAPACK.hetrs_rook!(B.uplo, B.LD, B.ipiv, R)
end
else
if issymmetric(B)
LAPACK.sytrs!(B.uplo, B.LD, B.ipiv, R)
else
LAPACK.hetrs!(B.uplo, B.LD, B.ipiv, R)
end
end
end
# There is no fallback solver for Bunch-Kaufman so we'll have to promote to same element type
function ldiv!(B::BunchKaufman{T}, R::StridedVecOrMat{S}) where {T,S}
TS = promote_type(T,S)
return ldiv!(convert(BunchKaufman{TS}, B), convert(AbstractArray{TS}, R))
end
function logabsdet(F::BunchKaufman)
M = F.LD
p = F.ipiv
n = size(F.LD, 1)
if !issuccess(F)
return eltype(F)(-Inf), zero(eltype(F))
end
s = one(real(eltype(F)))
i = 1
abs_det = zero(real(eltype(F)))
while i <= n
if p[i] > 0
elm = M[i,i]
s *= sign(elm)
abs_det += log(abs(elm))
i += 1
else
# 2x2 pivot case. Make sure not to square before the subtraction by scaling
# with the off-diagonal element. This is safe because the off diagonal is
# always large for 2x2 pivots.
if F.uplo == 'U'
elm = M[i, i + 1]*(M[i,i]/M[i, i + 1]*M[i + 1, i + 1] -
(issymmetric(F) ? M[i, i + 1] : conj(M[i, i + 1])))
s *= sign(elm)
abs_det += log(abs(elm))
else
elm = M[i + 1,i]*(M[i, i]/M[i + 1, i]*M[i + 1, i + 1] -
(issymmetric(F) ? M[i + 1, i] : conj(M[i + 1, i])))
s *= sign(elm)
abs_det += log(abs(elm))
end
i += 2
end
end
return abs_det, s
end
## reconstruct the original matrix
## TODO: understand the procedure described at
## http://www.nag.com/numeric/FL/nagdoc_fl22/pdf/F07/f07mdf.pdf