Skip to content

Rely on LinearAlgebra.eigsort! for sorting eigenvalues of the #151

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

Merged
merged 1 commit into from
Apr 20, 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
92 changes: 61 additions & 31 deletions src/eigenSelfAdjoint.jl
Original file line number Diff line number Diff line change
Expand Up @@ -162,7 +162,7 @@ function eig2x2!(d::StridedVector, e::StridedVector, j::Integer, vectors::Matrix
return c, s
end

function eigvalsPWK!(S::SymTridiagonal{T}; tol = eps(T), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) where {T<:Real}
function eigvalsPWK!(S::SymTridiagonal{T}; tol = eps(T)) where {T<:Real}
d = S.dv
e = S.ev
n = length(d)
Expand Down Expand Up @@ -217,9 +217,7 @@ function eigvalsPWK!(S::SymTridiagonal{T}; tol = eps(T), sortby::Union{Function,
end
end

# LinearAlgebra.eigvals will pass sortby=nothing but LAPACK always sort the symmetric
# eigenvalue problem so we'll will do the same here
LinearAlgebra.sorteig!(d, sortby === nothing ? LinearAlgebra.eigsortby : sortby)
return d
end

function eigQL!(
Expand Down Expand Up @@ -287,8 +285,8 @@ function eigQL!(
end
end
end
p = sortperm(d)
return d[p], vectors[:, p]

return d, vectors
end

function eigQR!(
Expand Down Expand Up @@ -339,8 +337,8 @@ function eigQR!(
end
end
end
p = sortperm(d)
return d[p], vectors[:, p]

return d, vectors
end

# Own implementation based on Parlett's book
Expand Down Expand Up @@ -575,45 +573,71 @@ function symtriUpper!(
)
end

_eigvals!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
LinearAlgebra.sorteig!(eigvalsPWK!(A; tol), sortby == nothing ? LinearAlgebra.eigsortby : sortby)

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

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

_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)
_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)

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

LinearAlgebra.eigvals!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigvals!(A; tol, sortby)
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)
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)...)
LinearAlgebra.Eigen(
LinearAlgebra.sorteig!(
eigQL!(
A.diagonals;
vectors = Array(A.Q),
tol
)...,
sortby == nothing ? LinearAlgebra.eigsortby : sortby
)...
)

_eigen!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = LinearAlgebra.Eigen(
eigQL!(A, vectors = Matrix{eltype(A)}(I, size(A, 1), size(A, 1)), tol = tol)...,
)
_eigen!(A::SymTridiagonal; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
LinearAlgebra.Eigen(
LinearAlgebra.sorteig!(
eigQL!(
A;
vectors = Matrix{eltype(A)}(I, size(A, 1), size(A, 1)),
tol
)...,
sortby == nothing ? LinearAlgebra.eigsortby : sortby
)...
)

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

_eigen!(A::Symmetric{<:Real}; tol = eps(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, sortby)

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

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

LinearAlgebra.eigen!(A::Hermitian; tol = eps(real(eltype(A))), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) = _eigen!(A; tol, sortby)
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)
LinearAlgebra.eigen!(A::Symmetric{<:Real}; tol = eps(eltype(A)), sortby::Union{Function,Nothing}=LinearAlgebra.eigsortby) =
_eigen!(A; tol, sortby)


function eigen2!(
Expand All @@ -623,29 +647,35 @@ function eigen2!(
V = zeros(eltype(A), 2, size(A, 1))
V[1] = 1
V[end] = 1
eigQL!(A.diagonals, vectors = rmul!(V, A.Q), tol = tol)
LinearAlgebra.sorteig!(
eigQL!(A.diagonals; vectors = rmul!(V, A.Q), tol)...,
LinearAlgebra.eigsortby
)
end

function eigen2!(A::SymTridiagonal; tol = eps(real(float(one(eltype(A))))))
V = zeros(eltype(A), 2, size(A, 1))
V[1] = 1
V[end] = 1
eigQL!(A, vectors = V, tol = tol)
LinearAlgebra.sorteig!(
eigQL!(A; vectors = V, tol)...,
LinearAlgebra.eigsortby
)
end

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

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


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

eigen2(A::Hermitian, 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)

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

# First method of each type here is identical to the method defined in
# LinearAlgebra but is needed for disambiguation
Expand Down
3 changes: 1 addition & 2 deletions test/eigenselfadjoint.jl
Original file line number Diff line number Diff line change
Expand Up @@ -26,9 +26,8 @@ using Test, GenericLinearAlgebra, LinearAlgebra, Quaternions
@testset "QR version (QL is default)" begin
vals, vecs =
GenericLinearAlgebra.eigQR!(copy(T), vectors = Matrix{eltype(T)}(I, n, n))
@test issorted(vals)
@test (vecs' * T) * vecs Diagonal(vals)
@test eigvals(T) vals
@test eigvals(T) sort(vals)
@test vecs'vecs Matrix(I, n, n)
end
end
Expand Down
Loading