-
-
Notifications
You must be signed in to change notification settings - Fork 5.5k
/
symmetriceigen.jl
281 lines (230 loc) · 10.2 KB
/
symmetriceigen.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
# This file is a part of Julia. License is MIT: https://julialang.org/license
# Eigensolvers for symmetric and Hermitian matrices
eigen!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}; sortby::Union{Function,Nothing}=nothing) =
Eigen(sorteig!(LAPACK.syevr!('V', 'A', A.uplo, A.data, 0.0, 0.0, 0, 0, -1.0)..., sortby)...)
function eigen(A::RealHermSymComplexHerm; sortby::Union{Function,Nothing}=nothing)
T = eltype(A)
S = eigtype(T)
eigen!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), sortby=sortby)
end
eigen!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}, irange::UnitRange) =
Eigen(LAPACK.syevr!('V', 'I', A.uplo, A.data, 0.0, 0.0, irange.start, irange.stop, -1.0)...)
"""
eigen(A::Union{SymTridiagonal, Hermitian, Symmetric}, irange::UnitRange) -> Eigen
Compute the eigenvalue decomposition of `A`, returning an [`Eigen`](@ref) factorization object `F`
which contains the eigenvalues in `F.values` and the eigenvectors in the columns of the
matrix `F.vectors`. (The `k`th eigenvector can be obtained from the slice `F.vectors[:, k]`.)
Iterating the decomposition produces the components `F.values` and `F.vectors`.
The following functions are available for `Eigen` objects: [`inv`](@ref), [`det`](@ref), and [`isposdef`](@ref).
The [`UnitRange`](@ref) `irange` specifies indices of the sorted eigenvalues to search for.
!!! note
If `irange` is not `1:n`, where `n` is the dimension of `A`, then the returned factorization
will be a *truncated* factorization.
"""
function eigen(A::RealHermSymComplexHerm, irange::UnitRange)
T = eltype(A)
S = eigtype(T)
eigen!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), irange)
end
eigen!(A::RealHermSymComplexHerm{T,<:StridedMatrix}, vl::Real, vh::Real) where {T<:BlasReal} =
Eigen(LAPACK.syevr!('V', 'V', A.uplo, A.data, convert(T, vl), convert(T, vh), 0, 0, -1.0)...)
"""
eigen(A::Union{SymTridiagonal, Hermitian, Symmetric}, vl::Real, vu::Real) -> Eigen
Compute the eigenvalue decomposition of `A`, returning an [`Eigen`](@ref) factorization object `F`
which contains the eigenvalues in `F.values` and the eigenvectors in the columns of the
matrix `F.vectors`. (The `k`th eigenvector can be obtained from the slice `F.vectors[:, k]`.)
Iterating the decomposition produces the components `F.values` and `F.vectors`.
The following functions are available for `Eigen` objects: [`inv`](@ref), [`det`](@ref), and [`isposdef`](@ref).
`vl` is the lower bound of the window of eigenvalues to search for, and `vu` is the upper bound.
!!! note
If [`vl`, `vu`] does not contain all eigenvalues of `A`, then the returned factorization
will be a *truncated* factorization.
"""
function eigen(A::RealHermSymComplexHerm, vl::Real, vh::Real)
T = eltype(A)
S = eigtype(T)
eigen!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), vl, vh)
end
function eigvals!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}; sortby::Union{Function,Nothing}=nothing)
vals = LAPACK.syevr!('N', 'A', A.uplo, A.data, 0.0, 0.0, 0, 0, -1.0)[1]
!isnothing(sortby) && sort!(vals, by=sortby)
return vals
end
function eigvals(A::RealHermSymComplexHerm; sortby::Union{Function,Nothing}=nothing)
T = eltype(A)
S = eigtype(T)
eigvals!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), sortby=sortby)
end
"""
eigvals!(A::Union{SymTridiagonal, Hermitian, Symmetric}, irange::UnitRange) -> values
Same as [`eigvals`](@ref), but saves space by overwriting the input `A`, instead of creating a copy.
`irange` is a range of eigenvalue *indices* to search for - for instance, the 2nd to 8th eigenvalues.
"""
eigvals!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}, irange::UnitRange) =
LAPACK.syevr!('N', 'I', A.uplo, A.data, 0.0, 0.0, irange.start, irange.stop, -1.0)[1]
"""
eigvals(A::Union{SymTridiagonal, Hermitian, Symmetric}, irange::UnitRange) -> values
Return the eigenvalues of `A`. It is possible to calculate only a subset of the
eigenvalues by specifying a [`UnitRange`](@ref) `irange` covering indices of the sorted eigenvalues,
e.g. the 2nd to 8th eigenvalues.
# Examples
```jldoctest
julia> A = SymTridiagonal([1.; 2.; 1.], [2.; 3.])
3×3 SymTridiagonal{Float64, Vector{Float64}}:
1.0 2.0 ⋅
2.0 2.0 3.0
⋅ 3.0 1.0
julia> eigvals(A, 2:2)
1-element Vector{Float64}:
0.9999999999999996
julia> eigvals(A)
3-element Vector{Float64}:
-2.1400549446402604
1.0000000000000002
5.140054944640259
```
"""
function eigvals(A::RealHermSymComplexHerm, irange::UnitRange)
T = eltype(A)
S = eigtype(T)
eigvals!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), irange)
end
"""
eigvals!(A::Union{SymTridiagonal, Hermitian, Symmetric}, vl::Real, vu::Real) -> values
Same as [`eigvals`](@ref), but saves space by overwriting the input `A`, instead of creating a copy.
`vl` is the lower bound of the interval to search for eigenvalues, and `vu` is the upper bound.
"""
eigvals!(A::RealHermSymComplexHerm{T,<:StridedMatrix}, vl::Real, vh::Real) where {T<:BlasReal} =
LAPACK.syevr!('N', 'V', A.uplo, A.data, convert(T, vl), convert(T, vh), 0, 0, -1.0)[1]
"""
eigvals(A::Union{SymTridiagonal, Hermitian, Symmetric}, vl::Real, vu::Real) -> values
Return the eigenvalues of `A`. It is possible to calculate only a subset of the eigenvalues
by specifying a pair `vl` and `vu` for the lower and upper boundaries of the eigenvalues.
# Examples
```jldoctest
julia> A = SymTridiagonal([1.; 2.; 1.], [2.; 3.])
3×3 SymTridiagonal{Float64, Vector{Float64}}:
1.0 2.0 ⋅
2.0 2.0 3.0
⋅ 3.0 1.0
julia> eigvals(A, -1, 2)
1-element Vector{Float64}:
1.0000000000000009
julia> eigvals(A)
3-element Vector{Float64}:
-2.1400549446402604
1.0000000000000002
5.140054944640259
```
"""
function eigvals(A::RealHermSymComplexHerm, vl::Real, vh::Real)
T = eltype(A)
S = eigtype(T)
eigvals!(S != T ? convert(AbstractMatrix{S}, A) : copy(A), vl, vh)
end
eigmax(A::RealHermSymComplexHerm{<:Real,<:StridedMatrix}) = eigvals(A, size(A, 1):size(A, 1))[1]
eigmin(A::RealHermSymComplexHerm{<:Real,<:StridedMatrix}) = eigvals(A, 1:1)[1]
function eigen!(A::HermOrSym{T,S}, B::HermOrSym{T,S}; sortby::Union{Function,Nothing}=nothing) where {T<:BlasReal,S<:StridedMatrix}
vals, vecs, _ = LAPACK.sygvd!(1, 'V', A.uplo, A.data, B.uplo == A.uplo ? B.data : copy(B.data'))
GeneralizedEigen(sorteig!(vals, vecs, sortby)...)
end
function eigen!(A::Hermitian{T,S}, B::Hermitian{T,S}; sortby::Union{Function,Nothing}=nothing) where {T<:BlasComplex,S<:StridedMatrix}
vals, vecs, _ = LAPACK.sygvd!(1, 'V', A.uplo, A.data, B.uplo == A.uplo ? B.data : copy(B.data'))
GeneralizedEigen(sorteig!(vals, vecs, sortby)...)
end
function eigen!(A::RealHermSymComplexHerm{T,S}, B::AbstractMatrix{T}; sortby::Union{Function,Nothing}=nothing) where {T<:Number,S<:StridedMatrix}
U = cholesky(B).U
vals, w = eigen!(UtiAUi!(A, U))
vecs = U \ w
GeneralizedEigen(sorteig!(vals, vecs, sortby)...)
end
# Perform U' \ A / U in-place.
UtiAUi!(As::Symmetric, Utr::UpperTriangular) = Symmetric(_UtiAsymUi!(As.uplo, parent(As), parent(Utr)), sym_uplo(As.uplo))
UtiAUi!(As::Hermitian, Utr::UpperTriangular) = Hermitian(_UtiAsymUi!(As.uplo, parent(As), parent(Utr)), sym_uplo(As.uplo))
UtiAUi!(As::Symmetric, Udi::Diagonal) = Symmetric(_UtiAsymUi_diag!(As.uplo, parent(As), Udi), sym_uplo(As.uplo))
UtiAUi!(As::Hermitian, Udi::Diagonal) = Hermitian(_UtiAsymUi_diag!(As.uplo, parent(As), Udi), sym_uplo(As.uplo))
# U is upper triangular
function _UtiAsymUi!(uplo, A, U)
n = size(A, 1)
μ⁻¹ = 1 / U[1, 1]
αμ⁻² = A[1, 1] * μ⁻¹' * μ⁻¹
# Update (1, 1) element
A[1, 1] = αμ⁻²
if n > 1
Unext = view(U, 2:n, 2:n)
if uplo === 'U'
# Update submatrix
for j in 2:n, i in 2:j
A[i, j] = (
A[i, j]
- μ⁻¹' * U[1, j] * A[1, i]'
- μ⁻¹ * A[1, j] * U[1, i]'
+ αμ⁻² * U[1, j] * U[1, i]'
)
end
# Update vector
for j in 2:n
A[1, j] = A[1, j] * μ⁻¹' - U[1, j] * αμ⁻²
end
ldiv!(view(A', 2:n, 1), UpperTriangular(Unext)', view(A', 2:n, 1))
else
# Update submatrix
for j in 2:n, i in 2:j
A[j, i] = (
A[j, i]
- μ⁻¹ * A[i, 1]' * U[1, j]'
- μ⁻¹' * U[1, i] * A[j, 1]
+ αμ⁻² * U[1, i] * U[1, j]'
)
end
# Update vector
for j in 2:n
A[j, 1] = A[j, 1] * μ⁻¹ - U[1, j]' * αμ⁻²
end
ldiv!(view(A, 2:n, 1), UpperTriangular(Unext)', view(A, 2:n, 1))
end
# Recurse
_UtiAsymUi!(uplo, view(A, 2:n, 2:n), Unext)
end
return A
end
# U is diagonal
function _UtiAsymUi_diag!(uplo, A, U)
n = size(A, 1)
μ⁻¹ = 1 / U[1, 1]
αμ⁻² = A[1, 1] * μ⁻¹' * μ⁻¹
# Update (1, 1) element
A[1, 1] = αμ⁻²
if n > 1
Unext = view(U, 2:n, 2:n)
if uplo === 'U'
# No need to update any submatrix when U is diagonal
# Update vector
for j in 2:n
A[1, j] = A[1, j] * μ⁻¹'
end
ldiv!(view(A', 2:n, 1), Diagonal(Unext)', view(A', 2:n, 1))
else
# No need to update any submatrix when U is diagonal
# Update vector
for j in 2:n
A[j, 1] = A[j, 1] * μ⁻¹
end
ldiv!(view(A, 2:n, 1), Diagonal(Unext)', view(A, 2:n, 1))
end
# Recurse
_UtiAsymUi!(uplo, view(A, 2:n, 2:n), Unext)
end
return A
end
function eigvals!(A::HermOrSym{T,S}, B::HermOrSym{T,S}; sortby::Union{Function,Nothing}=nothing) where {T<:BlasReal,S<:StridedMatrix}
vals = LAPACK.sygvd!(1, 'N', A.uplo, A.data, B.uplo == A.uplo ? B.data : copy(B.data'))[1]
isnothing(sortby) || sort!(vals, by=sortby)
return vals
end
function eigvals!(A::Hermitian{T,S}, B::Hermitian{T,S}; sortby::Union{Function,Nothing}=nothing) where {T<:BlasComplex,S<:StridedMatrix}
vals = LAPACK.sygvd!(1, 'N', A.uplo, A.data, B.uplo == A.uplo ? B.data : copy(B.data'))[1]
isnothing(sortby) || sort!(vals, by=sortby)
return vals
end
eigvecs(A::HermOrSym) = eigvecs(eigen(A))