From d0af9c528e5c551444bb31c866ff08f9e0481d78 Mon Sep 17 00:00:00 2001 From: Jishnu Bhattacharya Date: Fri, 8 Nov 2024 20:00:56 +0530 Subject: [PATCH] Avoid constprop in `syevd!` and `syev!` (#56442) This improves compilation times slightly: ```julia julia> using LinearAlgebra julia> A = rand(2,2); julia> @time eigen!(Hermitian(A)); 0.163380 seconds (180.51 k allocations: 8.760 MiB, 99.88% compilation time) # master 0.155285 seconds (163.77 k allocations: 7.971 MiB, 99.87% compilation time) # This PR ``` The idea is that the constant propagation is only required to infer the return type, and isn't necessary in the body of the method. We may therefore annotate the body with a `@constprop :none`. --- stdlib/LinearAlgebra/src/lapack.jl | 24 ++++++++++++++++-------- 1 file changed, 16 insertions(+), 8 deletions(-) diff --git a/stdlib/LinearAlgebra/src/lapack.jl b/stdlib/LinearAlgebra/src/lapack.jl index 97dff0031329b..5c2b66881585c 100644 --- a/stdlib/LinearAlgebra/src/lapack.jl +++ b/stdlib/LinearAlgebra/src/lapack.jl @@ -5318,6 +5318,14 @@ solution `X`. """ hetrs!(uplo::AbstractChar, A::AbstractMatrix, ipiv::AbstractVector{BlasInt}, B::AbstractVecOrMat) +for f in (:syevd!, :syev!) + _f = Symbol(:_, f) + @eval function $f(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix) + W, A = $_f(jobz, uplo, A) + jobz == 'V' ? (W, A) : W + end +end + # Symmetric (real) eigensolvers for (syev, syevr, syevd, sygvd, elty) in ((:dsyev_,:dsyevr_,:dsyevd_,:dsygvd_,:Float64), @@ -5329,7 +5337,7 @@ for (syev, syevr, syevd, sygvd, elty) in # INTEGER INFO, LDA, LWORK, N # * .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ) - Base.@constprop :aggressive function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + Base.@constprop :none function _syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) @chkvalidparam 1 jobz ('N', 'V') chkuplo(uplo) @@ -5350,7 +5358,7 @@ for (syev, syevr, syevd, sygvd, elty) in resize!(work, lwork) end end - jobz == 'V' ? (W, A) : W + W, A end # SUBROUTINE DSYEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, @@ -5429,7 +5437,7 @@ for (syev, syevr, syevd, sygvd, elty) in # * .. Array Arguments .. # INTEGER IWORK( * ) # DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ) - Base.@constprop :aggressive function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + Base.@constprop :none function _syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) @@ -5459,7 +5467,7 @@ for (syev, syevr, syevd, sygvd, elty) in resize!(iwork, liwork) end end - jobz == 'V' ? (W, A) : W + W, A end # Generalized eigenproblem @@ -5526,7 +5534,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # * .. Array Arguments .. # DOUBLE PRECISION RWORK( * ), W( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) - Base.@constprop :aggressive function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + Base.@constprop :none function _syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) @@ -5550,7 +5558,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in resize!(work, lwork) end end - jobz == 'V' ? (W, A) : W + W, A end # SUBROUTINE ZHEEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, @@ -5639,7 +5647,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # INTEGER IWORK( * ) # DOUBLE PRECISION RWORK( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) - Base.@constprop :aggressive function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + Base.@constprop :none function _syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) @@ -5673,7 +5681,7 @@ for (syev, syevr, syevd, sygvd, elty, relty) in resize!(iwork, liwork) end end - jobz == 'V' ? (W, A) : W + W, A end # SUBROUTINE ZHEGVD( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,