Skip to content

Commit 686804d

Browse files
authored
LAPACK: Aggressive constprop to concretely infer syev!/syevd! (#55295)
Currently, these are inferred as a 2-Tuple of possible return types depending on `jobz`, but since `jobz` is usually a constant, we may propagate it aggressively and have the return types inferred concretely.
1 parent 125bac4 commit 686804d

File tree

2 files changed

+14
-4
lines changed

2 files changed

+14
-4
lines changed

stdlib/LinearAlgebra/src/lapack.jl

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -5329,7 +5329,7 @@ for (syev, syevr, syevd, sygvd, elty) in
53295329
# INTEGER INFO, LDA, LWORK, N
53305330
# * .. Array Arguments ..
53315331
# DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
5332-
function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
5332+
Base.@constprop :aggressive function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
53335333
require_one_based_indexing(A)
53345334
@chkvalidparam 1 jobz ('N', 'V')
53355335
chkuplo(uplo)
@@ -5429,7 +5429,7 @@ for (syev, syevr, syevd, sygvd, elty) in
54295429
# * .. Array Arguments ..
54305430
# INTEGER IWORK( * )
54315431
# DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
5432-
function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
5432+
Base.@constprop :aggressive function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
54335433
require_one_based_indexing(A)
54345434
@chkvalidparam 1 jobz ('N', 'V')
54355435
chkstride1(A)
@@ -5526,7 +5526,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in
55265526
# * .. Array Arguments ..
55275527
# DOUBLE PRECISION RWORK( * ), W( * )
55285528
# COMPLEX*16 A( LDA, * ), WORK( * )
5529-
function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
5529+
Base.@constprop :aggressive function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
55305530
require_one_based_indexing(A)
55315531
@chkvalidparam 1 jobz ('N', 'V')
55325532
chkstride1(A)
@@ -5639,7 +5639,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in
56395639
# INTEGER IWORK( * )
56405640
# DOUBLE PRECISION RWORK( * )
56415641
# COMPLEX*16 A( LDA, * ), WORK( * )
5642-
function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
5642+
Base.@constprop :aggressive function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty})
56435643
require_one_based_indexing(A)
56445644
@chkvalidparam 1 jobz ('N', 'V')
56455645
chkstride1(A)

stdlib/LinearAlgebra/test/lapack.jl

Lines changed: 10 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -889,4 +889,14 @@ end
889889
@test UpperTriangular(A) == UpperTriangular(B)
890890
end
891891

892+
@testset "inference in syev!/syevd!" begin
893+
for T in (Float32, Float64), CT in (T, Complex{T})
894+
A = rand(CT, 4,4)
895+
@inferred (A -> LAPACK.syev!('N', 'U', A))(A)
896+
@inferred (A -> LAPACK.syev!('V', 'U', A))(A)
897+
@inferred (A -> LAPACK.syevd!('N', 'U', A))(A)
898+
@inferred (A -> LAPACK.syevd!('V', 'U', A))(A)
899+
end
900+
end
901+
892902
end # module TestLAPACK

0 commit comments

Comments
 (0)