diff --git a/base/exports.jl b/base/exports.jl index 09ea9466190cf..39148a90e55c9 100644 --- a/base/exports.jl +++ b/base/exports.jl @@ -566,6 +566,8 @@ export eig, eigfact, eigfact!, + eigmax, + eigmin, eigs, eigvals, eigvecs, diff --git a/base/linalg.jl b/base/linalg.jl index 974e4c1b33b6b..8c63cfe78be59 100644 --- a/base/linalg.jl +++ b/base/linalg.jl @@ -49,6 +49,8 @@ export eig, eigfact, eigfact!, + eigmax, + eigmin, eigs, eigvals, eigvecs, diff --git a/base/linalg/factorization.jl b/base/linalg/factorization.jl index 049e149a5b5d6..d1a177eeb0a19 100644 --- a/base/linalg/factorization.jl +++ b/base/linalg/factorization.jl @@ -443,6 +443,16 @@ end eigvals(x::Number) = [one(x)] +#Computes maximum and minimum eigenvalue +function eigmax(A::Union(Number, StridedMatrix)) + v = eigvals(A) + iscomplex(v) ? error("Complex eigenvalues cannot be ordered") : max(v) +end +function eigmin(A::Union(Number, StridedMatrix)) + v = eigvals(A) + iscomplex(v) ? error("Complex eigenvalues cannot be ordered") : min(v) +end + inv(A::EigenDense) = diagmm(A.vectors, 1.0/A.values)*A.vectors' det(A::EigenDense) = prod(A.values) diff --git a/base/linalg/hermitian.jl b/base/linalg/hermitian.jl index 97ab39d85acb6..0b2d4756f7db0 100644 --- a/base/linalg/hermitian.jl +++ b/base/linalg/hermitian.jl @@ -32,6 +32,7 @@ eigvals(A::Hermitian, il::Int, ih::Int) = LAPACK.syevr!('N', 'I', A.uplo, copy(A eigvals(A::Hermitian, vl::Real, vh::Real) = LAPACK.syevr!('N', 'V', A.uplo, copy(A.S), vl, vh, 0, 0, -1.0)[1] eigvals(A::Hermitian) = eigvals(A, 1, size(A, 1)) eigmax(A::Hermitian) = eigvals(A, size(A, 1), size(A, 1))[1] +eigmin(A::Hermitian) = eigvals(A, 1, 1)[1] function expm(A::Hermitian) F = eigfact(A) diff --git a/base/linalg/tridiag.jl b/base/linalg/tridiag.jl index c79ee7e588e0a..a26a1884e6791 100644 --- a/base/linalg/tridiag.jl +++ b/base/linalg/tridiag.jl @@ -67,6 +67,11 @@ eig(m::SymTridiagonal) = LAPACK.stegr!('V', copy(m.dv), copy(m.ev)) eigvals(m::SymTridiagonal, il::Int, iu::Int) = LAPACK.stebz!('I', 'E', 0.0, 0.0, il, iu, -1.0, copy(m.dv), copy(m.ev))[1] eigvals(m::SymTridiagonal, vl::Float64, vu::Float64) = LAPACK.stebz!('V', 'E', vl, vu, 0, 0, -1.0, copy(m.dv), copy(m.ev))[1] eigvals(m::SymTridiagonal) = LAPACK.stebz!('A', 'E', 0.0, 0.0, 0, 0, -1.0, copy(m.dv), copy(m.ev))[1] + +#Computes largest and smallest eigenvalue +eigmax(m::SymTridiagonal) = eigvals(m, size(m, 1), size(m, 1))[1] +eigmin(m::SymTridiagonal) = eigvals(m, 1, 1)[1] + #Compute selected eigenvectors only corresponding to particular eigenvalues eigvecs(m::SymTridiagonal) = eig(m)[2] eigvecs{Eigenvalue<:Real}(m::SymTridiagonal, eigvals::Vector{Eigenvalue}) = LAPACK.stein!(m.dv, m.ev, eigvals) diff --git a/doc/stdlib/linalg.rst b/doc/stdlib/linalg.rst index bc9799ddb7285..241eb614c7c28 100644 --- a/doc/stdlib/linalg.rst +++ b/doc/stdlib/linalg.rst @@ -95,11 +95,19 @@ Linear algebra functions in Julia are largely implemented by calling functions f Returns the eigenvalues of ``A``. +.. function:: eigmax(A) + + Returns the largest eigenvalue of ``A``. + +.. function:: eigmin(A) + + Returns the smallest eigenvalue of ``A``. + .. function:: eigvecs(A, [eigvals]) Returns the eigenvectors of ``A``. - If the optional vector of eigenvalues ``eigvals`` is specified, returns the specific corresponding eigenvectors. (Currently this optional syntax only works for SymTridiagonal matrices.) + For SymTridiagonal matrices, if the optional vector of eigenvalues ``eigvals`` is specified, returns the specific corresponding eigenvectors. .. function:: eigfact(A)