From ed3fd1eeaa7d7db5ca6a0eac037c994cd79712a3 Mon Sep 17 00:00:00 2001 From: Aravindh Krishnamoorthy Date: Tue, 5 Dec 2023 20:21:08 +0100 Subject: [PATCH] Update the bunchkaufman.jl/getproperties! function to align to the anti-aliasing implementation in #51763 + include comments in symmetriceigen.jl --- stdlib/LinearAlgebra/src/bunchkaufman.jl | 6 +++--- stdlib/LinearAlgebra/src/symmetriceigen.jl | 3 +++ 2 files changed, 6 insertions(+), 3 deletions(-) diff --git a/stdlib/LinearAlgebra/src/bunchkaufman.jl b/stdlib/LinearAlgebra/src/bunchkaufman.jl index fa464815ac5c4..f995b3b8444c4 100644 --- a/stdlib/LinearAlgebra/src/bunchkaufman.jl +++ b/stdlib/LinearAlgebra/src/bunchkaufman.jl @@ -279,7 +279,7 @@ Base.propertynames(B::BunchKaufman, private::Bool=false) = (:p, :P, :L, :U, :D, (private ? fieldnames(typeof(B)) : ())...) function getproperties!(B::BunchKaufman{T,<:StridedMatrix}) where {T<:BlasFloat} - # NOTE: Copy of BunchKaufman's getproperty(B, :D) and getproperty(B, :L/U) since 'getproperty' is not in place. + # NOTE: Unlike in the 'getproperty' function, in this function L/U and D are computed in place. if B.rook LUD, od = LAPACK.syconvf_rook!(B.uplo, 'C', B.LD, B.ipiv) else @@ -289,12 +289,12 @@ function getproperties!(B::BunchKaufman{T,<:StridedMatrix}) where {T<:BlasFloat} M = UnitUpperTriangular(LUD) du = od[2:end] # Avoid aliasing dl and du. - dl = B.symmetric ? copy(du) : conj.(du) + dl = B.symmetric ? du : conj.(du) else M = UnitLowerTriangular(LUD) dl = od[1:end-1] # Avoid aliasing dl and du. - du = B.symmetric ? copy(dl) : conj.(dl) + du = B.symmetric ? dl : conj.(dl) end return (M, Tridiagonal(dl, diag(LUD), du), B.p) end diff --git a/stdlib/LinearAlgebra/src/symmetriceigen.jl b/stdlib/LinearAlgebra/src/symmetriceigen.jl index f78092a2780ad..4e35616f62181 100644 --- a/stdlib/LinearAlgebra/src/symmetriceigen.jl +++ b/stdlib/LinearAlgebra/src/symmetriceigen.jl @@ -286,6 +286,9 @@ end function eigvals!(A::StridedMatrix{T}, F::LU{T,<:StridedMatrix}; sortby::Union{Function,Nothing}=nothing) where {T} L = UnitLowerTriangular(F.L) U = UpperTriangular(F.U) + # Compute generalized eigenvalues of equivalent matrix: + # A' = inv(L)*(P*A)*inv(U) + # See: https://github.com/JuliaLang/julia/pull/50471#issuecomment-1627836781 permuterows!(A, F.p) ldiv!(L, A) rdiv!(A, U)