Skip to content

Commit 41db513

Browse files
authored
Eigvecs for specific eigvals for Symmetric/Hermitian (#1268)
This implements the suggestion in #1248. After this, ```julia julia> S = Symmetric(rand(3,3)) 3×3 Symmetric{Float64, Matrix{Float64}}: 0.376244 0.895193 0.332219 0.895193 0.563134 0.148036 0.332219 0.148036 0.711689 julia> vals = eigvals(S) 3-element Vector{Float64}: -0.45018177966363415 0.5911683834592292 1.5100812026842658 julia> eigvecs(S, vals[1:2]) 3×2 Matrix{Float64}: 0.752343 -0.16258 -0.645224 -0.37737 -0.132912 0.911679 ```
1 parent c33045c commit 41db513

File tree

3 files changed

+27
-2
lines changed

3 files changed

+27
-2
lines changed

src/symmetriceigen.jl

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -331,6 +331,12 @@ function eigvals!(A::Hermitian{T,S}, B::Hermitian{T,S}; sortby::Union{Function,N
331331
end
332332
eigvecs(A::HermOrSym) = eigvecs(eigen(A))
333333

334+
function eigvecs(A::RealHermSymComplexHerm, eigvals::AbstractVector{<:Real})
335+
F = hessenberg(A) # transform to SymTridiagonal form
336+
X = eigvecs(F.H, eigvals)
337+
return F.Q * X # transform eigvecs of F.H back to eigvecs of A
338+
end
339+
334340
function eigvals(A::AbstractMatrix, C::Cholesky; sortby::Union{Function,Nothing}=nothing)
335341
if ishermitian(A)
336342
eigvals!(eigencopy_oftype(Hermitian(A), eigtype(eltype(A))), C; sortby)

src/tridiag.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -310,7 +310,7 @@ eigmin(A::SymTridiagonal) = eigvals(A, 1:1)[1]
310310

311311
#Compute selected eigenvectors only corresponding to particular eigenvalues
312312
"""
313-
eigvecs(A::SymTridiagonal[, eigvals]) -> Matrix
313+
eigvecs(A::Union{Symmetric, Hermitian, SymTridiagonal}[, eigvals])::Matrix
314314
315315
Return a matrix `M` whose columns are the eigenvectors of `A`. (The `k`th eigenvector can
316316
be obtained from the slice `M[:, k]`.)
@@ -345,7 +345,7 @@ julia> eigvecs(A, [1.])
345345
-0.5547001962252291
346346
```
347347
"""
348-
eigvecs(A::SymTridiagonal{<:BlasFloat,<:StridedVector}, eigvals::Vector{<:Real}) = LAPACK.stein!(A.dv, A.ev, eigvals)
348+
eigvecs(A::SymTridiagonal{<:BlasFloat,<:StridedVector}, eigvals::StridedVector{<:Real}) = LAPACK.stein!(A.dv, A.ev, eigvals)
349349

350350
function svdvals!(A::SymTridiagonal)
351351
vals = eigvals!(A)

test/symmetriceigen.jl

Lines changed: 19 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -201,4 +201,23 @@ end
201201
end
202202
end
203203

204+
@testset "eigvecs for specific eigvals" begin
205+
function testeigvecs(S, vals)
206+
V = eigvecs(S, vals)
207+
@test S * V V * Diagonal(vals)
208+
end
209+
for T in (Symmetric, Hermitian)
210+
S = T(rand(3,3))
211+
vals = eigvals(S)
212+
testeigvecs(S, vals)
213+
testeigvecs(S, vals[1:2])
214+
testeigvecs(S, @view vals[2:2])
215+
end
216+
H = Hermitian(rand(ComplexF64,3,3))
217+
vals = eigvals(H)
218+
testeigvecs(H, vals)
219+
testeigvecs(H, vals[1:2])
220+
testeigvecs(H, @view vals[2:2])
221+
end
222+
204223
end # module TestSymmetricEigen

0 commit comments

Comments
 (0)