Skip to content

Commit dce8c1c

Browse files
Improve accuracy of 2x2 symmetric/Hermitian eigendecomposition (#1158)
* Improve accuracy of 2x2 symmetric/Hermitian eigendecomposition * bump version * Update src/eigen.jl Co-authored-by: Yuto Horikawa <[email protected]> --------- Co-authored-by: Yuto Horikawa <[email protected]>
1 parent b4e49ce commit dce8c1c

File tree

3 files changed

+8
-10
lines changed

3 files changed

+8
-10
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "StaticArrays"
22
uuid = "90137ffa-7385-5640-81b9-e52037218182"
3-
version = "1.5.23"
3+
version = "1.5.24"
44

55
[deps]
66
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"

src/eigen.jl

Lines changed: 2 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -159,14 +159,10 @@ end
159159
@inline function _eig(::Size{(2,2)}, A::LinearAlgebra.RealHermSymComplexHerm{T}, permute, scale) where {T <: Real}
160160
a = A.data
161161
TA = eltype(A)
162-
163162
@inbounds if A.uplo == 'U'
164163
if !iszero(a[3]) # A is not diagonal
165164
t_half = real(a[1] + a[4]) / 2
166-
d = real(a[1] * a[4] - a[3]' * a[3]) # Should be real
167-
168-
tmp2 = t_half * t_half - d
169-
tmp = tmp2 < 0 ? zero(tmp2) : sqrt(tmp2) # Numerically stable for identity matrices, etc.
165+
tmp = norm(SVector((a[1] - a[4])/2, a[3]'))
170166
vals = SVector(t_half - tmp, t_half + tmp)
171167

172168
v11 = vals[1] - a[4]
@@ -187,10 +183,7 @@ end
187183
else # A.uplo == 'L'
188184
if !iszero(a[2]) # A is not diagonal
189185
t_half = real(a[1] + a[4]) / 2
190-
d = real(a[1] * a[4] - a[2]' * a[2]) # Should be real
191-
192-
tmp2 = t_half * t_half - d
193-
tmp = tmp2 < 0 ? zero(tmp2) : sqrt(tmp2) # Numerically stable for identity matrices, etc.
186+
tmp = norm(SVector((a[1] - a[4])/2, a[2]))
194187
vals = SVector(t_half - tmp, t_half + tmp)
195188

196189
v11 = vals[1] - a[4]

test/eigen.jl

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -93,6 +93,11 @@ using StaticArrays, Test, LinearAlgebra
9393
@test (@inferred SVector{2,ComplexF64} eigvals(Symmetric(m1), Symmetric(m2))) eigvals(Symmetric(m1_a), Symmetric(m2_a))
9494
end
9595

96+
@testset "2×2, difficult case" begin
97+
m1 = SA[1.6590891025248637 -2.7087777909606814e-7; -2.7087777909606814e-7 1.659089317183428]
98+
@test norm(m1 - Array(eigen(m1))) < 1e-15
99+
end
100+
96101
@test_throws DimensionMismatch eigvals(SA[1 2 3; 4 5 6], SA[1 2 3; 4 5 5])
97102
@test_throws DimensionMismatch eigvals(SA[1 2; 4 5], SA[1 2 3; 4 5 5; 3 4 5])
98103

0 commit comments

Comments
 (0)