Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

add methods for symmetric real eigen #145

Merged
merged 2 commits into from
Jan 15, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 11 additions & 2 deletions src/eigenSelfAdjoint.jl
Original file line number Diff line number Diff line change
Expand Up @@ -448,8 +448,7 @@ function singleShiftQR!(
end

symtri!(A::Hermitian) = A.uplo == 'L' ? symtriLower!(A.data) : symtriUpper!(A.data)
symtri!(A::Symmetric{T}) where {T<:Real} =
A.uplo == 'L' ? symtriLower!(A.data) : symtriUpper!(A.data)
symtri!(A::Symmetric{<:Real}) = A.uplo == 'L' ? symtriLower!(A.data) : symtriUpper!(A.data)

# Assume that lower triangle stores the relevant part
function symtriLower!(
Expand Down Expand Up @@ -580,6 +579,7 @@ _eigvals!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Union{Function,

_eigvals!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = eigvals!(symtri!(A); tol, sortby)

_eigvals!(A::Symmetric{<:Real}; tol = eps(eltype(A)), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = eigvals!(symtri!(A); tol, sortby)

LinearAlgebra.eigvals!(A::SymmetricTridiagonalFactorization; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
_eigvals!(A; tol, sortby)
Expand All @@ -589,6 +589,7 @@ LinearAlgebra.eigvals!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Un

LinearAlgebra.eigvals!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigvals!(A; tol, sortby)

LinearAlgebra.eigvals!(A::Symmetric{<:Real}; tol = eps(eltype(A)), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigvals!(A; tol, sortby)

_eigen!(A::SymmetricTridiagonalFactorization; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
LinearAlgebra.Eigen(LinearAlgebra.sorteig!(eigQL!(A.diagonals, vectors = Array(A.Q), tol = tol)..., sortby)...)
Expand All @@ -599,6 +600,7 @@ _eigen!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Union{Function,No

_eigen!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigen!(symtri!(A), tol = tol)

_eigen!(A::Symmetric{<:Real}; tol = eps(eltype(A)), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigen!(symtri!(A), tol = tol)

LinearAlgebra.eigen!(A::SymmetricTridiagonalFactorization; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
_eigen!(A; tol, sortby)
Expand All @@ -607,6 +609,8 @@ LinearAlgebra.eigen!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Unio

LinearAlgebra.eigen!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigen!(A; tol, sortby)

LinearAlgebra.eigen!(A::Symmetric{<:Real}; tol = eps(eltype(A)), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigen!(A; tol, sortby)


function eigen2!(
A::SymmetricTridiagonalFactorization;
Expand All @@ -628,12 +632,17 @@ end
eigen2!(A::Hermitian; tol = eps(float(real(one(eltype(A)))))) =
eigen2!(symtri!(A), tol = tol)

eigen2!(A::Symmetric{<:Real}; tol = eps(float(one(eltype(A))))) =
eigen2!(symtri!(A), tol = tol)


eigen2(A::SymTridiagonal; tol = eps(float(real(one(eltype(A)))))) =
eigen2!(copy(A), tol = tol)

eigen2(A::Hermitian, tol = eps(float(real(one(eltype(A)))))) = eigen2!(copy(A), tol = tol)

eigen2(A::Symmetric{<:Real}, tol = eps(float(one(eltype(A))))) = eigen2!(copy(A), tol = tol)

# First method of each type here is identical to the method defined in
# LinearAlgebra but is needed for disambiguation
const _eigencopy_oftype = if VERSION >= v"1.9"
Expand Down
2 changes: 1 addition & 1 deletion test/eigenselfadjoint.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ using Test, GenericLinearAlgebra, LinearAlgebra, Quaternions
end

@testset "(full) Symmetric" for uplo in (:L, :U)
A = Hermitian(big.(randn(n, n)), uplo)
A = Symmetric(big.(randn(n, n)), uplo)
vals, vecs = eigen(A)
@testset "default" begin
@test vecs' * A * vecs ≈ diagm(0 => vals)
Expand Down
Loading