From dcd1fb2148dd7affee9449e304a9b42e7feba1ce Mon Sep 17 00:00:00 2001 From: Jishnu Bhattacharya Date: Tue, 12 Mar 2024 07:34:25 +0530 Subject: [PATCH] LAPACK: validate input parameters to throw informative errors (#53631) This PR validates the input parameters to the Julia LAPACK wrappers, so that the error messages are more informative. On nightly ```julia julia> using LinearAlgebra julia> LAPACK.geev!('X', 'X', rand(2,2)) ** On entry to DGEEV parameter number 1 had an illegal value ERROR: ArgumentError: invalid argument #1 to LAPACK call ``` This PR ```julia julia> using LinearAlgebra julia> LAPACK.geev!('X', 'X', rand(2,2)) ERROR: ArgumentError: argument #1: jobvl must be one of ('N', 'V'), but 'X' was passed ``` Secondly, moved certain allocations (e.g. in `geevx`) below the validation checks, so that these only happen for valid parameter values. Thirdly, added `require_one_based_indexing` checks to functions where these were missing. --- stdlib/LinearAlgebra/src/blas.jl | 2 +- stdlib/LinearAlgebra/src/lapack.jl | 240 ++++++++++++++++++++++------ stdlib/LinearAlgebra/test/lapack.jl | 184 ++++++++++++++++----- 3 files changed, 344 insertions(+), 82 deletions(-) diff --git a/stdlib/LinearAlgebra/src/blas.jl b/stdlib/LinearAlgebra/src/blas.jl index 7bf7f8bf95100..8903aaddcc11b 100644 --- a/stdlib/LinearAlgebra/src/blas.jl +++ b/stdlib/LinearAlgebra/src/blas.jl @@ -168,7 +168,7 @@ end "Check that upper/lower (for special matrices) is correctly specified" function chkuplo(uplo::AbstractChar) if !(uplo == 'U' || uplo == 'L') - throw(ArgumentError(lazy"uplo argument must be 'U' (upper) or 'L' (lower), got $uplo")) + throw(ArgumentError(lazy"uplo argument must be 'U' (upper) or 'L' (lower), got '$uplo'")) end uplo end diff --git a/stdlib/LinearAlgebra/src/lapack.jl b/stdlib/LinearAlgebra/src/lapack.jl index a85978435bf46..02dfa0079038b 100644 --- a/stdlib/LinearAlgebra/src/lapack.jl +++ b/stdlib/LinearAlgebra/src/lapack.jl @@ -55,10 +55,23 @@ function chkposdef(ret::BlasInt) end end +# Generic fallback function to assert that parameters are valid +# In specific cases, the following functions may be more useful +macro chkvalidparam(position::Int, param, validvalues) + :(chkvalidparam($position, $(string(param)), $(esc(param)), $validvalues)) +end +function chkvalidparam(position::Int, var::String, val, validvals) + if val ∉ validvals + throw(ArgumentError( + "argument #$position: $var must be one of $validvals, but $(repr(val)) was passed")) + end + return val +end + "Check that {c}transpose is correctly specified" function chktrans(trans::AbstractChar) if !(trans == 'N' || trans == 'C' || trans == 'T') - throw(ArgumentError("trans argument must be 'N' (no transpose), 'T' (transpose), or 'C' (conjugate transpose), got $trans")) + throw(ArgumentError("trans argument must be 'N' (no transpose), 'T' (transpose), or 'C' (conjugate transpose), got '$trans'")) end trans end @@ -66,7 +79,7 @@ end "Check that left/right hand side multiply is correctly specified" function chkside(side::AbstractChar) if !(side == 'L' || side == 'R') - throw(ArgumentError("side argument must be 'L' (left hand multiply) or 'R' (right hand multiply), got $side")) + throw(ArgumentError("side argument must be 'L' (left hand multiply) or 'R' (right hand multiply), got '$side'")) end side end @@ -74,7 +87,7 @@ end "Check that unit diagonal flag is correctly specified" function chkdiag(diag::AbstractChar) if !(diag == 'U' || diag =='N') - throw(ArgumentError("diag argument must be 'U' (unit diagonal) or 'N' (non-unit diagonal), got $diag")) + throw(ArgumentError("diag argument must be 'U' (unit diagonal) or 'N' (non-unit diagonal), got '$diag'")) end diag end @@ -93,6 +106,7 @@ end function chkuplofinite(A::AbstractMatrix, uplo::AbstractChar) require_one_based_indexing(A) + chkuplo(uplo) m, n = size(A) if uplo == 'U' @inbounds for j in 1:n, i in 1:j @@ -213,7 +227,9 @@ for (gebal, gebak, elty, relty) in # .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ), SCALE( * ) function gebal!(job::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) + @chkvalidparam 1 job ('N', 'P', 'S', 'B') n = checksquare(A) chkfinite(A) # balancing routines don't support NaNs and Infs ihi = Ref{BlasInt}() @@ -238,6 +254,7 @@ for (gebal, gebak, elty, relty) in ilo::BlasInt, ihi::BlasInt, scale::AbstractVector{$relty}, V::AbstractMatrix{$elty}) require_one_based_indexing(scale, V) + @chkvalidparam 1 job ('N', 'P', 'S', 'B') chkstride1(scale, V) chkside(side) chkfinite(V) # balancing routines don't support NaNs and Infs @@ -624,11 +641,13 @@ does not equal zero then the `j`th column of `A` is permuted to the front of geqp3!(A::AbstractMatrix, jpvt::AbstractVector{BlasInt}, tau::AbstractVector) function geqp3!(A::AbstractMatrix{<:BlasFloat}, jpvt::AbstractVector{BlasInt}) + require_one_based_indexing(A, jpvt) m, n = size(A) geqp3!(A, jpvt, similar(A, min(m, n))) end function geqp3!(A::AbstractMatrix{<:BlasFloat}) + require_one_based_indexing(A) m, n = size(A) geqp3!(A, zeros(BlasInt, n), similar(A, min(m, n))) end @@ -777,6 +796,7 @@ for (larfg, elty) in # .. Array Arguments .. # DOUBLE PRECISION x( * ) function larfg!(x::AbstractVector{$elty}) + require_one_based_indexing(x) N = BlasInt(length(x)) α = Ref{$elty}(x[1]) incx = BlasInt(1) @@ -805,6 +825,7 @@ for (larf, elty) in # DOUBLE PRECISION c( ldc, * ), v( * ), work( * ) function larf!(side::AbstractChar, v::AbstractVector{$elty}, τ::$elty, C::AbstractMatrix{$elty}, work::AbstractVector{$elty}) + require_one_based_indexing(v, C, work) m, n = size(C) chkside(side) ldc = max(1, stride(C, 2)) @@ -820,6 +841,7 @@ for (larf, elty) in function larf!(side::AbstractChar, v::AbstractVector{$elty}, τ::$elty, C::AbstractMatrix{$elty}) + require_one_based_indexing(v, C) m, n = size(C) chkside(side) lwork = side == 'L' ? n : m @@ -1134,6 +1156,7 @@ for (gesvx, elty) in AF::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}, equed::AbstractChar, R::AbstractVector{$elty}, C::AbstractVector{$elty}, B::AbstractVecOrMat{$elty}) require_one_based_indexing(A, AF, ipiv, R, C, B) + @chkvalidparam 1 fact ('F', 'N', 'E') chktrans(trans) chkstride1(ipiv, R, C, B) n = checksquare(A) @@ -1168,6 +1191,7 @@ for (gesvx, elty) in end function gesvx!(A::AbstractMatrix{$elty}, B::AbstractVecOrMat{$elty}) + require_one_based_indexing(A, B) n = size(A,1) X, equed, R, C, B, rcond, ferr, berr, rpgf = gesvx!('N', 'N', A, @@ -1204,6 +1228,7 @@ for (gesvx, elty, relty) in AF::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}, equed::AbstractChar, R::AbstractVector{$relty}, C::AbstractVector{$relty}, B::AbstractVecOrMat{$elty}) require_one_based_indexing(A, AF, ipiv, R, C, B) + @chkvalidparam 1 fact ('F', 'N', 'E') chktrans(trans) chkstride1(A, AF, ipiv, R, C, B) n = checksquare(A) @@ -1239,6 +1264,7 @@ for (gesvx, elty, relty) in #Wrapper for the no-equilibration, no-transpose calculation function gesvx!(A::AbstractMatrix{$elty}, B::AbstractVecOrMat{$elty}) + require_one_based_indexing(A, B) n = size(A,1) X, equed, R, C, B, rcond, ferr, berr, rpgf = gesvx!('N', 'N', A, @@ -1577,8 +1603,11 @@ for (geev, gesvd, gesdd, ggsvd, elty, relty) in # DOUBLE PRECISION A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ), # $ WI( * ), WORK( * ), WR( * ) function geev!(jobvl::AbstractChar, jobvr::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) + @chkvalidparam 1 jobvl ('N', 'V') + @chkvalidparam 2 jobvr ('N', 'V') chkfinite(A) # balancing routines don't support NaNs and Infs lvecs = jobvl == 'V' rvecs = jobvr == 'V' @@ -1635,6 +1664,7 @@ for (geev, gesvd, gesdd, ggsvd, elty, relty) in function gesdd!(job::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) chkstride1(A) + @chkvalidparam 1 job ('A', 'S', 'O', 'N') m, n = size(A) minmn = min(m, n) if job == 'A' @@ -1716,6 +1746,9 @@ for (geev, gesvd, gesdd, ggsvd, elty, relty) in function gesvd!(jobu::AbstractChar, jobvt::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) chkstride1(A) + @chkvalidparam 1 jobu ('A', 'S', 'O', 'N') + @chkvalidparam 2 jobvt ('A', 'S', 'O', 'N') + (jobu == jobvt == 'O') && throw(ArgumentError("jobu and jobvt cannot both be O")) m, n = size(A) minmn = min(m, n) S = similar(A, $relty, minmn) @@ -1785,6 +1818,9 @@ for (geev, gesvd, gesdd, ggsvd, elty, relty) in function ggsvd!(jobu::AbstractChar, jobv::AbstractChar, jobq::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) require_one_based_indexing(A, B) chkstride1(A, B) + @chkvalidparam 1 jobu ('U', 'N') + @chkvalidparam 2 jobv ('V', 'N') + @chkvalidparam 3 jobq ('Q', 'N') m, n = size(A) if size(B, 2) != n throw(DimensionMismatch("B has second dimension $(size(B,2)) but needs $n")) @@ -1912,6 +1948,9 @@ for (f, elty) in ((:dggsvd3_, :Float64), function ggsvd3!(jobu::AbstractChar, jobv::AbstractChar, jobq::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) require_one_based_indexing(A, B) chkstride1(A, B) + @chkvalidparam 1 jobu ('U', 'N') + @chkvalidparam 2 jobv ('V', 'N') + @chkvalidparam 3 jobq ('Q', 'N') m, n = size(A) if size(B, 2) != n throw(DimensionMismatch("B has second dimension $(size(B,2)) but needs $n")) @@ -1971,6 +2010,9 @@ for (f, elty, relty) in ((:zggsvd3_, :ComplexF64, :Float64), function ggsvd3!(jobu::AbstractChar, jobv::AbstractChar, jobq::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) require_one_based_indexing(A, B) chkstride1(A, B) + @chkvalidparam 1 jobu ('U', 'N') + @chkvalidparam 2 jobv ('V', 'N') + @chkvalidparam 3 jobq ('Q', 'N') m, n = size(A) if size(B, 2) != n throw(DimensionMismatch("B has second dimension $(size(B,2)) but needs $n")) @@ -2057,14 +2099,13 @@ for (geevx, ggev, ggev3, elty) in # $ SCALE( * ), VL( LDVL, * ), VR( LDVR, * ), # $ WI( * ), WORK( * ), WR( * ) function geevx!(balanc::AbstractChar, jobvl::AbstractChar, jobvr::AbstractChar, sense::AbstractChar, A::AbstractMatrix{$elty}) - n = checksquare(A) - chkfinite(A) # balancing routines don't support NaNs and Infs - lda = max(1,stride(A,2)) - wr = similar(A, $elty, n) - wi = similar(A, $elty, n) - if balanc ∉ ['N', 'P', 'S', 'B'] - throw(ArgumentError("balanc must be 'N', 'P', 'S', or 'B', but $balanc was passed")) + require_one_based_indexing(A) + @chkvalidparam 1 balanc ('N', 'P', 'S', 'B') + @chkvalidparam 4 sense ('N', 'E', 'V', 'B') + if sense ∈ ('E', 'B') && !(jobvl == jobvr == 'V') + throw(ArgumentError("sense = '$sense' requires jobvl = 'V' and jobvr = 'V'")) end + n = checksquare(A) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2073,7 +2114,6 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - VL = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2082,6 +2122,11 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end + chkfinite(A) # balancing routines don't support NaNs and Infs + lda = max(1,stride(A,2)) + wr = similar(A, $elty, n) + wi = similar(A, $elty, n) + VL = similar(A, $elty, ldvl, n) VR = similar(A, $elty, ldvr, n) ilo = Ref{BlasInt}() ihi = Ref{BlasInt}() @@ -2143,11 +2188,6 @@ for (geevx, ggev, ggev3, elty) in if n != m throw(DimensionMismatch("A has dimensions $(size(A)), and B has dimensions $(size(B)), but A and B must have the same size")) end - lda = max(1, stride(A, 2)) - ldb = max(1, stride(B, 2)) - alphar = similar(A, $elty, n) - alphai = similar(A, $elty, n) - beta = similar(A, $elty, n) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2156,7 +2196,6 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - vl = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2165,6 +2204,12 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end + lda = max(1, stride(A, 2)) + ldb = max(1, stride(B, 2)) + alphar = similar(A, $elty, n) + alphai = similar(A, $elty, n) + beta = similar(A, $elty, n) + vl = similar(A, $elty, ldvl, n) vr = similar(A, $elty, ldvr, n) work = Vector{$elty}(undef, 1) lwork = BlasInt(-1) @@ -2207,11 +2252,6 @@ for (geevx, ggev, ggev3, elty) in if n != m throw(DimensionMismatch("A has dimensions $(size(A)), and B has dimensions $(size(B)), but A and B must have the same size")) end - lda = max(1, stride(A, 2)) - ldb = max(1, stride(B, 2)) - alphar = similar(A, $elty, n) - alphai = similar(A, $elty, n) - beta = similar(A, $elty, n) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2220,7 +2260,6 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - vl = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2229,6 +2268,12 @@ for (geevx, ggev, ggev3, elty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end + lda = max(1, stride(A, 2)) + ldb = max(1, stride(B, 2)) + alphar = similar(A, $elty, n) + alphai = similar(A, $elty, n) + beta = similar(A, $elty, n) + vl = similar(A, $elty, ldvl, n) vr = similar(A, $elty, ldvr, n) work = Vector{$elty}(undef, 1) lwork = BlasInt(-1) @@ -2275,13 +2320,17 @@ for (geevx, ggev, ggev3, elty, relty) in # COMPLEX*16 A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ), # $ W( * ), WORK( * ) function geevx!(balanc::AbstractChar, jobvl::AbstractChar, jobvr::AbstractChar, sense::AbstractChar, A::AbstractMatrix{$elty}) - n = checksquare(A) - chkfinite(A) # balancing routines don't support NaNs and Infs - lda = max(1,stride(A,2)) - w = similar(A, $elty, n) - if balanc ∉ ['N', 'P', 'S', 'B'] + require_one_based_indexing(A) + if balanc ∉ ('N', 'P', 'S', 'B') throw(ArgumentError("balanc must be 'N', 'P', 'S', or 'B', but $balanc was passed")) end + if sense ∉ ('N','E','V','B') + throw(ArgumentError("sense must be 'N', 'E', 'V' or 'B', but $sense was passed")) + end + if sense ∈ ('E', 'B') && !(jobvl == jobvr == 'V') + throw(ArgumentError("sense = '$sense' requires jobvl = 'V' and jobvr = 'V'")) + end + n = checksquare(A) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2290,7 +2339,6 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - VL = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2299,9 +2347,10 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end - if sense ∉ ['N','E','V','B'] - throw(ArgumentError("sense must be 'N', 'E', 'V' or 'B', but $sense was passed")) - end + chkfinite(A) # balancing routines don't support NaNs and Infs + lda = max(1,stride(A,2)) + w = similar(A, $elty, n) + VL = similar(A, $elty, ldvl, n) VR = similar(A, $elty, ldvr, n) ilo = Ref{BlasInt}() ihi = Ref{BlasInt}() @@ -2354,10 +2403,6 @@ for (geevx, ggev, ggev3, elty, relty) in if n != m throw(DimensionMismatch("A has dimensions $(size(A)), and B has dimensions $(size(B)), but A and B must have the same size")) end - lda = max(1, stride(A, 2)) - ldb = max(1, stride(B, 2)) - alpha = similar(A, $elty, n) - beta = similar(A, $elty, n) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2366,7 +2411,6 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - vl = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2375,6 +2419,11 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end + lda = max(1, stride(A, 2)) + ldb = max(1, stride(B, 2)) + alpha = similar(A, $elty, n) + beta = similar(A, $elty, n) + vl = similar(A, $elty, ldvl, n) vr = similar(A, $elty, ldvr, n) work = Vector{$elty}(undef, 1) lwork = BlasInt(-1) @@ -2419,10 +2468,6 @@ for (geevx, ggev, ggev3, elty, relty) in if n != m throw(DimensionMismatch("A has dimensions $(size(A)), and B has dimensions $(size(B)), but A and B must have the same size")) end - lda = max(1, stride(A, 2)) - ldb = max(1, stride(B, 2)) - alpha = similar(A, $elty, n) - beta = similar(A, $elty, n) ldvl = 0 if jobvl == 'V' ldvl = n @@ -2431,7 +2476,6 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvl must be 'V' or 'N', but $jobvl was passed")) end - vl = similar(A, $elty, ldvl, n) ldvr = 0 if jobvr == 'V' ldvr = n @@ -2440,6 +2484,11 @@ for (geevx, ggev, ggev3, elty, relty) in else throw(ArgumentError("jobvr must be 'V' or 'N', but $jobvr was passed")) end + lda = max(1, stride(A, 2)) + ldb = max(1, stride(B, 2)) + alpha = similar(A, $elty, n) + beta = similar(A, $elty, n) + vl = similar(A, $elty, ldvl, n) vr = similar(A, $elty, ldvr, n) work = Vector{$elty}(undef, 1) lwork = BlasInt(-1) @@ -2524,6 +2573,7 @@ for (laic1, elty) in function laic1!(job::Integer, x::AbstractVector{$elty}, sest::$elty, w::AbstractVector{$elty}, gamma::$elty) require_one_based_indexing(x, w) + @chkvalidparam 1 job (1,2) j = length(x) if j != length(w) throw(DimensionMismatch("vectors must have same length, but length of x is $j and length of w is $(length(w))")) @@ -2558,6 +2608,7 @@ for (laic1, elty, relty) in function laic1!(job::Integer, x::AbstractVector{$elty}, sest::$relty, w::AbstractVector{$elty}, gamma::$elty) require_one_based_indexing(x, w) + @chkvalidparam 1 job (1,2) j = length(x) if j != length(w) throw(DimensionMismatch("vectors must have same length, but length of x is $j and length of w is $(length(w))")) @@ -3303,6 +3354,7 @@ for (posv, potrf, potri, potrs, pstrf, elty, rtyp) in # DOUBLE PRECISION A( LDA, * ), WORK( 2*N ) # INTEGER PIV( N ) function pstrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}, tol::Real) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -3534,6 +3586,7 @@ for (trtri, trtrs, elty) in # .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ) function trtri!(uplo::AbstractChar, diag::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -3616,10 +3669,12 @@ for (trcon, trevc, trrfs, elty) in # INTEGER IWORK( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function trcon!(norm::AbstractChar, uplo::AbstractChar, diag::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) chkdiag(diag) n = checksquare(A) chkuplo(uplo) + @chkvalidparam 1 norm ('O', '1', 'I') rcond = Ref{$elty}() work = Vector{$elty}(undef, 3n) iwork = Vector{BlasInt}(undef, n) @@ -3651,9 +3706,10 @@ for (trcon, trevc, trrfs, elty) in VR::AbstractMatrix{$elty} = similar(T)) require_one_based_indexing(select, T, VL, VR) # Extract - if side ∉ ['L','R','B'] + if side ∉ ('L','R','B') throw(ArgumentError("side argument must be 'L' (left eigenvectors), 'R' (right eigenvectors), or 'B' (both), got $side")) end + @chkvalidparam 2 howmny ('A', 'B', 'S') n, mm = checksquare(T), size(VL, 2) ldt, ldvl, ldvr = stride(T, 2), stride(VL, 2), stride(VR, 2) @@ -3749,8 +3805,10 @@ for (trcon, trevc, trrfs, elty, relty) in # DOUBLE PRECISION RWORK( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function trcon!(norm::AbstractChar, uplo::AbstractChar, diag::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) + @chkvalidparam 1 norm ('O', '1', 'I') chkuplo(uplo) chkdiag(diag) rcond = Ref{$relty}(1) @@ -3790,9 +3848,10 @@ for (trcon, trevc, trrfs, elty, relty) in # Check chkstride1(T, select, VL, VR) - if side ∉ ['L','R','B'] + if side ∉ ('L','R','B') throw(ArgumentError("side argument must be 'L' (left eigenvectors), 'R' (right eigenvectors), or 'B' (both), got $side")) end + @chkvalidparam 2 howmny ('A', 'B', 'S') # Allocate m = Ref{BlasInt}() @@ -3919,6 +3978,7 @@ for (stev, stebz, stegr, stein, elty) in @eval begin function stev!(job::AbstractChar, dv::AbstractVector{$elty}, ev::AbstractVector{$elty}) require_one_based_indexing(dv, ev) + @chkvalidparam 1 job ('N', 'V') chkstride1(dv, ev) n = length(dv) if length(ev) != n - 1 && length(ev) != n @@ -3941,6 +4001,8 @@ for (stev, stebz, stegr, stein, elty) in #* eigenvalues. function stebz!(range::AbstractChar, order::AbstractChar, vl::$elty, vu::$elty, il::Integer, iu::Integer, abstol::Real, dv::AbstractVector{$elty}, ev::AbstractVector{$elty}) require_one_based_indexing(dv, ev) + @chkvalidparam 1 range ('A', 'V', 'I') + @chkvalidparam 2 order ('B', 'E') chkstride1(dv, ev) n = length(dv) if length(ev) != n - 1 @@ -3972,6 +4034,8 @@ for (stev, stebz, stegr, stein, elty) in function stegr!(jobz::AbstractChar, range::AbstractChar, dv::AbstractVector{$elty}, ev::AbstractVector{$elty}, vl::Real, vu::Real, il::Integer, iu::Integer) require_one_based_indexing(dv, ev) + @chkvalidparam 1 jobz ('N', 'V') + @chkvalidparam 2 range ('A', 'V', 'I') chkstride1(dv, ev) n = length(dv) ne = length(ev) @@ -4144,6 +4208,7 @@ for (syconv, sysv, sytrf, sytri, sytrs, elty) in # INTEGER IPIV( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function syconv!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -4201,6 +4266,7 @@ for (syconv, sysv, sytrf, sytri, sytrs, elty) in # INTEGER IPIV( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function sytrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -4225,6 +4291,8 @@ for (syconv, sysv, sytrf, sytri, sytrs, elty) in end function sytrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + chkuplo(uplo) n = checksquare(A) ipiv = similar(A, BlasInt, n) sytrf!(uplo, A, ipiv) @@ -4267,6 +4335,7 @@ for (syconv, sysv, sytrf, sytri, sytrs, elty) in # INTEGER IPIV( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function sytri!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -4358,6 +4427,7 @@ for (sysv, sytrf, sytri, sytrs, syconvf, elty) in # INTEGER IPIV( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function sytrf_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -4390,6 +4460,7 @@ for (sysv, sytrf, sytri, sytrs, syconvf, elty) in # INTEGER IPIV( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function sytri_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -4492,6 +4563,7 @@ for (syconv, hesv, hetrf, hetri, hetrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function syconv!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A,ipiv) chkstride1(A,ipiv) n = checksquare(A) chkuplo(uplo) @@ -4549,6 +4621,7 @@ for (syconv, hesv, hetrf, hetri, hetrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function hetrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -4570,6 +4643,8 @@ for (syconv, hesv, hetrf, hetri, hetrs, elty, relty) in end function hetrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + chkuplo(uplo) n = checksquare(A) ipiv = similar(A, BlasInt, n) hetrf!(uplo, A, ipiv) @@ -4614,6 +4689,7 @@ for (syconv, hesv, hetrf, hetri, hetrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function hetri!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -4638,6 +4714,7 @@ for (syconv, hesv, hetrf, hetri, hetrs, elty, relty) in function hetrs!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}, B::AbstractVecOrMat{$elty}) require_one_based_indexing(A, ipiv, B) + chkuplo(uplo) chkstride1(A,B,ipiv) n = checksquare(A) if n != size(B,1) @@ -4702,6 +4779,7 @@ for (hesv, hetrf, hetri, hetrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function hetrf_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -4732,6 +4810,7 @@ for (hesv, hetrf, hetri, hetrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function hetri_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A,ipiv) chkstride1(A,ipiv) n = checksquare(A) chkuplo(uplo) @@ -4757,6 +4836,7 @@ for (hesv, hetrf, hetri, hetrs, elty, relty) in ipiv::AbstractVector{BlasInt}, B::AbstractVecOrMat{$elty}) require_one_based_indexing(A, ipiv, B) chkstride1(A,B,ipiv) + chkuplo(uplo) n = checksquare(A) if n != size(B,1) throw(DimensionMismatch("B has first dimension $(size(B,1)), but needs $n")) @@ -4822,6 +4902,7 @@ for (sysv, sytrf, sytri, sytrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function sytrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -4846,6 +4927,8 @@ for (sysv, sytrf, sytri, sytrs, elty, relty) in end function sytrf!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + chkuplo(uplo) n = checksquare(A) ipiv = similar(A, BlasInt, n) sytrf!(uplo, A, ipiv) @@ -4889,6 +4972,7 @@ for (sysv, sytrf, sytri, sytrs, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function sytri!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -4980,6 +5064,7 @@ for (sysv, sytrf, sytri, sytrs, syconvf, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function sytrf_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -5013,6 +5098,7 @@ for (sysv, sytrf, sytri, sytrs, syconvf, elty, relty) in # INTEGER IPIV( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function sytri_rook!(uplo::AbstractChar, A::AbstractMatrix{$elty}, ipiv::AbstractVector{BlasInt}) + require_one_based_indexing(A, ipiv) chkstride1(A, ipiv) n = checksquare(A) chkuplo(uplo) @@ -5240,6 +5326,9 @@ for (syev, syevr, syevd, sygvd, elty) in # * .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ) function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') + chkuplo(uplo) chkstride1(A) n = checksquare(A) W = similar(A, $elty, n) @@ -5273,15 +5362,18 @@ for (syev, syevr, syevd, sygvd, elty) in # DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * ) function syevr!(jobz::AbstractChar, range::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}, vl::AbstractFloat, vu::AbstractFloat, il::Integer, iu::Integer, abstol::AbstractFloat) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') + @chkvalidparam 2 range ('A', 'V', 'I') chkstride1(A) n = checksquare(A) - chkuplofinite(A, uplo) if range == 'I' && !(1 <= il <= iu <= n) throw(ArgumentError("illegal choice of eigenvalue indices (il = $il, iu = $iu), which must be between 1 and n = $n")) end if range == 'V' && vl >= vu throw(ArgumentError("lower boundary, $vl, must be less than upper boundary, $vu")) end + chkuplofinite(A, uplo) lda = stride(A,2) m = Ref{BlasInt}() W = similar(A, $elty, n) @@ -5334,6 +5426,8 @@ for (syev, syevr, syevd, sygvd, elty) in # INTEGER IWORK( * ) # DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ) function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) n = checksquare(A) chkuplofinite(A, uplo) @@ -5375,6 +5469,10 @@ for (syev, syevr, syevd, sygvd, elty) in # INTEGER IWORK( * ) # DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * ) function sygvd!(itype::Integer, jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 itype 1:3 + @chkvalidparam 2 jobz ('N', 'V') + chkuplo(uplo) chkstride1(A, B) n, m = checksquare(A, B) if n != m @@ -5425,6 +5523,8 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # DOUBLE PRECISION RWORK( * ), W( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function syev!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) chkuplofinite(A, uplo) n = checksquare(A) @@ -5464,6 +5564,9 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # COMPLEX*16 A( LDA, * ), WORK( * ), Z( LDZ, * ) function syevr!(jobz::AbstractChar, range::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}, vl::AbstractFloat, vu::AbstractFloat, il::Integer, iu::Integer, abstol::AbstractFloat) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') + @chkvalidparam 2 range ('A', 'V', 'I') chkstride1(A) chkuplofinite(A, uplo) n = checksquare(A) @@ -5533,6 +5636,8 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # DOUBLE PRECISION RWORK( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function syevd!(jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) + @chkvalidparam 1 jobz ('N', 'V') chkstride1(A) chkuplofinite(A, uplo) n = checksquare(A) @@ -5578,6 +5683,9 @@ for (syev, syevr, syevd, sygvd, elty, relty) in # DOUBLE PRECISION RWORK( * ), W( * ) # COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * ) function sygvd!(itype::Integer, jobz::AbstractChar, uplo::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 itype 1:3 + @chkvalidparam 2 jobz ('N', 'V') chkstride1(A, B) chkuplofinite(A, uplo) chkuplofinite(B, uplo) @@ -5827,6 +5935,8 @@ for (gecon, elty) in # INTEGER IWORK( * ) # DOUBLE PRECISION A( LDA, * ), WORK( * ) function gecon!(normtype::AbstractChar, A::AbstractMatrix{$elty}, anorm::$elty) + require_one_based_indexing(A) + @chkvalidparam 1 normtype ('0', '1', 'I') chkstride1(A) n = checksquare(A) lda = max(1, stride(A, 2)) @@ -5861,6 +5971,8 @@ for (gecon, elty, relty) in # DOUBLE PRECISION RWORK( * ) # COMPLEX*16 A( LDA, * ), WORK( * ) function gecon!(normtype::AbstractChar, A::AbstractMatrix{$elty}, anorm::$relty) + require_one_based_indexing(A) + @chkvalidparam 1 normtype ('0', '1', 'I') chkstride1(A) n = checksquare(A) lda = max(1, stride(A, 2)) @@ -5903,6 +6015,7 @@ for (gehrd, elty) in # * .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * ) function gehrd!(ilo::Integer, ihi::Integer, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkfinite(A) # balancing routines don't support NaNs and Infs @@ -6005,6 +6118,8 @@ for (ormhr, elty) in require_one_based_indexing(A, tau, C) chkstride1(A, tau, C) + chkside(side) + chktrans(trans) n = checksquare(A) mC, nC = size(C, 1), size(C, 2) @@ -6052,6 +6167,8 @@ for (hseqr, elty) in function hseqr!(job::AbstractChar, compz::AbstractChar, ilo::Integer, ihi::Integer, H::AbstractMatrix{$elty}, Z::AbstractMatrix{$elty}) require_one_based_indexing(H, Z) + @chkvalidparam 1 job ('E', 'S') + @chkvalidparam 2 compz ('N', 'I', 'V') chkstride1(H) n = checksquare(H) checksquare(Z) == n || throw(DimensionMismatch()) @@ -6094,6 +6211,8 @@ for (hseqr, elty) in function hseqr!(job::AbstractChar, compz::AbstractChar, ilo::Integer, ihi::Integer, H::AbstractMatrix{$elty}, Z::AbstractMatrix{$elty}) require_one_based_indexing(H, Z) + @chkvalidparam 1 job ('E', 'S') + @chkvalidparam 2 compz ('N', 'I', 'V') chkstride1(H) n = checksquare(H) checksquare(Z) == n || throw(DimensionMismatch()) @@ -6152,6 +6271,7 @@ for (hetrd, elty) in # * .. Array Arguments .. # DOUBLE PRECISION A( LDA, * ), D( * ), E( * ), TAU( * ), WORK( * ) function hetrd!(uplo::AbstractChar, A::AbstractMatrix{$elty}) + require_one_based_indexing(A) chkstride1(A) n = checksquare(A) chkuplo(uplo) @@ -6257,7 +6377,9 @@ for (ormtr, elty) in require_one_based_indexing(A, tau, C) chkstride1(A, tau, C) n = checksquare(A) + chkside(side) chkuplo(uplo) + chktrans(trans) mC, nC = size(C, 1), size(C, 2) if n - length(tau) != 1 @@ -6305,6 +6427,7 @@ for (gees, gges, gges3, elty) in # $ WR( * ) function gees!(jobvs::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) + @chkvalidparam 1 jobvs ('N', 'V') chkstride1(A) n = checksquare(A) sdim = Vector{BlasInt}(undef, 1) @@ -6344,6 +6467,9 @@ for (gees, gges, gges3, elty) in # $ B( LDB, * ), BETA( * ), VSL( LDVSL, * ), # $ VSR( LDVSR, * ), WORK( * ) function gges!(jobvsl::AbstractChar, jobvsr::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 jobvsl ('N', 'V') + @chkvalidparam 2 jobvsr ('N', 'V') chkstride1(A, B) n, m = checksquare(A, B) if n != m @@ -6393,6 +6519,9 @@ for (gees, gges, gges3, elty) in # $ B( LDB, * ), BETA( * ), VSL( LDVSL, * ), # $ VSR( LDVSR, * ), WORK( * ) function gges3!(jobvsl::AbstractChar, jobvsr::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 jobvsl ('N', 'V') + @chkvalidparam 2 jobvsr ('N', 'V') chkstride1(A, B) n, m = checksquare(A, B) if n != m @@ -6448,6 +6577,7 @@ for (gees, gges, gges3, elty, relty) in # COMPLEX*16 A( LDA, * ), VS( LDVS, * ), W( * ), WORK( * ) function gees!(jobvs::AbstractChar, A::AbstractMatrix{$elty}) require_one_based_indexing(A) + @chkvalidparam 1 jobvs ('N', 'V') chkstride1(A) n = checksquare(A) sort = 'N' @@ -6489,6 +6619,9 @@ for (gees, gges, gges3, elty, relty) in # $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ), # $ WORK( * ) function gges!(jobvsl::AbstractChar, jobvsr::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 jobvsl ('N', 'V') + @chkvalidparam 2 jobvsr ('N', 'V') chkstride1(A, B) n, m = checksquare(A, B) if n != m @@ -6539,6 +6672,9 @@ for (gees, gges, gges3, elty, relty) in # $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ), # $ WORK( * ) function gges3!(jobvsl::AbstractChar, jobvsr::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}) + require_one_based_indexing(A, B) + @chkvalidparam 1 jobvsl ('N', 'V') + @chkvalidparam 2 jobvsr ('N', 'V') chkstride1(A, B) n, m = checksquare(A, B) if n != m @@ -6627,6 +6763,8 @@ for (trexc, trsen, tgsen, elty) in # * .. Array Arguments .. # DOUBLE PRECISION Q( LDQ, * ), T( LDT, * ), WORK( * ) function trexc!(compq::AbstractChar, ifst::BlasInt, ilst::BlasInt, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}) + require_one_based_indexing(T, Q) + @chkvalidparam 1 compq ('V', 'N') chkstride1(T, Q) n = checksquare(T) ldt = max(1, stride(T, 2)) @@ -6659,6 +6797,9 @@ for (trexc, trsen, tgsen, elty) in # DOUBLE PRECISION Q( LDQ, * ), T( LDT, * ), WI( * ), WORK( * ), WR( * ) function trsen!(job::AbstractChar, compq::AbstractChar, select::AbstractVector{BlasInt}, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}) + require_one_based_indexing(T, Q, select) + @chkvalidparam 1 job ('N', 'E', 'V', 'B') + @chkvalidparam 2 compq ('V', 'N') chkstride1(T, Q, select) n = checksquare(T) ldt = max(1, stride(T, 2)) @@ -6714,6 +6855,7 @@ for (trexc, trsen, tgsen, elty) in # .. function tgsen!(select::AbstractVector{BlasInt}, S::AbstractMatrix{$elty}, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}, Z::AbstractMatrix{$elty}) + require_one_based_indexing(select, S, T, Q, Z) chkstride1(select, S, T, Q, Z) n, nt, nq, nz = checksquare(S, T, Q, Z) if n != nt @@ -6779,6 +6921,8 @@ for (trexc, trsen, tgsen, elty, relty) in # .. Array Arguments .. # DOUBLE PRECISION Q( LDQ, * ), T( LDT, * ), WORK( * ) function trexc!(compq::AbstractChar, ifst::BlasInt, ilst::BlasInt, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}) + require_one_based_indexing(T, Q) + @chkvalidparam 1 compq ('V', 'N') chkstride1(T, Q) n = checksquare(T) ldt = max(1, stride(T, 2)) @@ -6809,6 +6953,9 @@ for (trexc, trsen, tgsen, elty, relty) in # COMPLEX Q( LDQ, * ), T( LDT, * ), W( * ), WORK( * ) function trsen!(job::AbstractChar, compq::AbstractChar, select::AbstractVector{BlasInt}, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}) + require_one_based_indexing(select, T, Q) + @chkvalidparam 1 job ('N', 'E', 'V', 'B') + @chkvalidparam 2 compq ('N', 'V') chkstride1(select, T, Q) n = checksquare(T) ldt = max(1, stride(T, 2)) @@ -6859,6 +7006,7 @@ for (trexc, trsen, tgsen, elty, relty) in # .. function tgsen!(select::AbstractVector{BlasInt}, S::AbstractMatrix{$elty}, T::AbstractMatrix{$elty}, Q::AbstractMatrix{$elty}, Z::AbstractMatrix{$elty}) + require_one_based_indexing(select, S, T, Q, Z) chkstride1(select, S, T, Q, Z) n, nt, nq, nz = checksquare(S, T, Q, Z) if n != nt @@ -6959,6 +7107,8 @@ for (fn, elty, relty) in ((:dtrsyl_, :Float64, :Float64), function trsyl!(transa::AbstractChar, transb::AbstractChar, A::AbstractMatrix{$elty}, B::AbstractMatrix{$elty}, C::AbstractMatrix{$elty}, isgn::Int=1) require_one_based_indexing(A, B, C) + chktrans(transa) + chktrans(transb) chkstride1(A, B, C) m, n = checksquare(A), checksquare(B) lda = max(1, stride(A, 2)) diff --git a/stdlib/LinearAlgebra/test/lapack.jl b/stdlib/LinearAlgebra/test/lapack.jl index 6e12c85204a78..652c6c2e27e6c 100644 --- a/stdlib/LinearAlgebra/test/lapack.jl +++ b/stdlib/LinearAlgebra/test/lapack.jl @@ -9,6 +9,7 @@ using LinearAlgebra: BlasInt @test_throws ArgumentError LinearAlgebra.LAPACK.chkside('Z') @test_throws ArgumentError LinearAlgebra.LAPACK.chkdiag('Z') @test_throws ArgumentError LinearAlgebra.LAPACK.chktrans('Z') +@test_throws ArgumentError LinearAlgebra.LAPACK.chkvalidparam(1, "job", 2, (0,1)) @testset "syevr" begin Random.seed!(123) @@ -35,6 +36,12 @@ using LinearAlgebra: BlasInt @test vals_test ≈ vals @test Z_test*(Diagonal(vals)*Z_test') ≈ Asym @test_throws DimensionMismatch LAPACK.sygvd!(1, 'V', 'U', copy(Asym), zeros(elty, 6, 6)) + + @test_throws "jobz must be one of ('N', 'V'), but 'X' was passed" LAPACK.syevr!('X', Asym) + @test_throws "jobz must be one of ('N', 'V'), but 'X' was passed" LAPACK.syev!('X', 'U', Asym) + @test_throws "uplo argument must be 'U' (upper) or 'L' (lower), got 'M'" LAPACK.syev!('N', 'M', Asym) + @test_throws "jobz must be one of ('N', 'V'), but 'X' was passed" LAPACK.syevd!('X', 'U', Asym) + @test_throws "uplo argument must be 'U' (upper) or 'L' (lower), got 'M'" LAPACK.syevd!('N', 'M', Asym) end end @@ -112,7 +119,9 @@ end D = LAPACK.gbtrs!('N',2,1,6,AB,ipiv,D) A = diagm(-2 => dl2, -1 => dl, 0 => d, 1 => du) @test A\C ≈ D - @test_throws DimensionMismatch LAPACK.gbtrs!('N',2,1,6,AB,ipiv,Matrix{elty}(undef,7,6)) + M = Matrix{elty}(undef,7,6) + @test_throws DimensionMismatch LAPACK.gbtrs!('N',2,1,6,AB,ipiv,M) + @test_throws ArgumentError LAPACK.gbtrs!('M',2,1,6,AB,ipiv,M) @test_throws LinearAlgebra.LAPACKException LAPACK.gbtrf!(2,1,6,zeros(elty,6,6)) end end @@ -141,9 +150,11 @@ end x10, x11 = Vector{LinearAlgebra.BlasInt}.(undef, (10, 11)) @test_throws DimensionMismatch LAPACK.gels!('N',A10x10,B11x11) @test_throws DimensionMismatch LAPACK.gels!('T',A10x10,B11x11) + @test_throws ArgumentError LAPACK.gels!('X',A10x10,B11x11) @test_throws DimensionMismatch LAPACK.gesv!(A10x10,B11x11) @test_throws DimensionMismatch LAPACK.getrs!('N',A10x10,x10,B11x11) @test_throws DimensionMismatch LAPACK.getrs!('T',A10x10,x10,B11x11) + @test_throws ArgumentError LAPACK.getrs!('X',A10x10,x10,B11x11) @test_throws DimensionMismatch LAPACK.getri!(A10x10,x11) end end @@ -177,12 +188,20 @@ end @test U ≈ lU @test S ≈ lS @test V' ≈ lVt + @test_throws ArgumentError LAPACK.gesvd!('X','S',A) + @test_throws ArgumentError LAPACK.gesvd!('S','X',A) B = rand(elty,10,10) # xggsvd3 replaced xggsvd in LAPACK 3.6.0 if LAPACK.version() < v"3.6.0" - @test_throws DimensionMismatch LAPACK.ggsvd!('S','S','S',A,B) + @test_throws DimensionMismatch LAPACK.ggsvd!('N','N','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd!('X','N','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd!('N','X','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd!('N','N','X',A,B) else - @test_throws DimensionMismatch LAPACK.ggsvd3!('S','S','S',A,B) + @test_throws DimensionMismatch LAPACK.ggsvd3!('N','N','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd3!('X','N','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd3!('N','X','N',A,B) + @test_throws ArgumentError LAPACK.ggsvd3!('N','N','X',A,B) end end end @@ -224,6 +243,7 @@ end X = rand(elty,10) B,Y,z = LAPACK.gels!('N',copy(A),copy(X)) @test A\X ≈ Y + @test_throws ArgumentError LAPACK.gels!('X',A,X) end end @@ -252,6 +272,9 @@ end fA = eigen(A, sortby=nothing) @test fA.values ≈ Aw @test fA.vectors ≈ Avr + + @test_throws ArgumentError LAPACK.geev!('X','V',A) + @test_throws ArgumentError LAPACK.geev!('N','X',A) end end @@ -284,6 +307,7 @@ end @test_throws DimensionMismatch LAPACK.gttrs!('N', x11, d, du, x9, y10, b) @test_throws DimensionMismatch LAPACK.gttrs!('N', dl, d, x11, x9, y10, b) @test_throws DimensionMismatch LAPACK.gttrs!('N', dl, d, du, x9, y10, x11) + @test_throws ArgumentError LAPACK.gttrs!('X', dl, d, du, x9, y10, x11) A = lu(Tridiagonal(dl,d,du)) b = rand(elty,10,5) c = copy(b) @@ -298,10 +322,17 @@ end A = rand(elty,10,10) A,tau = LAPACK.gelqf!(A) @test_throws DimensionMismatch LAPACK.orglq!(A,tau,11) - @test_throws DimensionMismatch LAPACK.ormlq!('R','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormlq!('L','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormlq!('R','N',A,zeros(elty,11),rand(elty,10,10)) - @test_throws DimensionMismatch LAPACK.ormlq!('L','N',A,zeros(elty,11),rand(elty,10,10)) + temp = rand(elty,11,11) + @test_throws DimensionMismatch LAPACK.ormlq!('R','N',A,tau,temp) + @test_throws DimensionMismatch LAPACK.ormlq!('L','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormlq!('X','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormlq!('R','X',A,tau,temp) + temp = zeros(elty,11) + B = copy(A) + @test_throws DimensionMismatch LAPACK.ormlq!('R','N',A,temp,B) + @test_throws DimensionMismatch LAPACK.ormlq!('L','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormlq!('X','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormlq!('L','X',A,temp,B) B = copy(A) C = LAPACK.orglq!(B,tau) @@ -312,30 +343,51 @@ end @test_throws DimensionMismatch LAPACK.orgqr!(A,tau,11) B = copy(A) @test LAPACK.orgqr!(B,tau) ≈ LAPACK.ormqr!('R','N',A,tau,Matrix{elty}(I, 10, 10)) - @test_throws DimensionMismatch LAPACK.ormqr!('R','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormqr!('L','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormqr!('R','N',A,zeros(elty,11),rand(elty,10,10)) - @test_throws DimensionMismatch LAPACK.ormqr!('L','N',A,zeros(elty,11),rand(elty,10,10)) + temp = rand(elty,11,11) + @test_throws DimensionMismatch LAPACK.ormqr!('R','N',A,tau,temp) + @test_throws DimensionMismatch LAPACK.ormqr!('L','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormqr!('X','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormqr!('L','X',A,tau,temp) + B = copy(A) + temp = zeros(elty,11) + @test_throws DimensionMismatch LAPACK.ormqr!('R','N',A,temp,B) + @test_throws DimensionMismatch LAPACK.ormqr!('L','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormqr!('X','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormqr!('L','X',A,temp,B) A = rand(elty,10,10) A,tau = LAPACK.geqlf!(A) @test_throws DimensionMismatch LAPACK.orgql!(A,tau,11) B = copy(A) @test LAPACK.orgql!(B,tau) ≈ LAPACK.ormql!('R','N',A,tau,Matrix{elty}(I, 10, 10)) - @test_throws DimensionMismatch LAPACK.ormql!('R','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormql!('L','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormql!('R','N',A,zeros(elty,11),rand(elty,10,10)) - @test_throws DimensionMismatch LAPACK.ormql!('L','N',A,zeros(elty,11),rand(elty,10,10)) + temp = rand(elty,11,11) + @test_throws DimensionMismatch LAPACK.ormql!('R','N',A,tau,temp) + @test_throws DimensionMismatch LAPACK.ormql!('L','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormql!('X','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormql!('L','X',A,tau,temp) + temp = zeros(elty,11) + B = copy(A) + @test_throws DimensionMismatch LAPACK.ormql!('R','N',A,temp,B) + @test_throws DimensionMismatch LAPACK.ormql!('L','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormql!('X','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormql!('L','X',A,temp,B) A = rand(elty,10,10) A,tau = LAPACK.gerqf!(A) @test_throws DimensionMismatch LAPACK.orgrq!(A,tau,11) B = copy(A) @test LAPACK.orgrq!(B,tau) ≈ LAPACK.ormrq!('R','N',A,tau,Matrix{elty}(I, 10, 10)) - @test_throws DimensionMismatch LAPACK.ormrq!('R','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormrq!('L','N',A,tau,rand(elty,11,11)) - @test_throws DimensionMismatch LAPACK.ormrq!('R','N',A,zeros(elty,11),rand(elty,10,10)) - @test_throws DimensionMismatch LAPACK.ormrq!('L','N',A,zeros(elty,11),rand(elty,10,10)) + temp = rand(elty,11,11) + @test_throws DimensionMismatch LAPACK.ormrq!('R','N',A,tau,temp) + @test_throws DimensionMismatch LAPACK.ormrq!('L','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormrq!('X','N',A,tau,temp) + @test_throws ArgumentError LAPACK.ormrq!('L','X',A,tau,temp) + B = copy(A) + temp = zeros(elty,11) + @test_throws DimensionMismatch LAPACK.ormrq!('R','N',A,temp,B) + @test_throws DimensionMismatch LAPACK.ormrq!('L','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormrq!('X','N',A,temp,B) + @test_throws ArgumentError LAPACK.ormrq!('L','X',A,temp,B) A = rand(elty,10,11) Q = copy(A) @@ -351,21 +403,29 @@ end T = zeros(elty,10,11) @test_throws DimensionMismatch LAPACK.gemqrt!('L','N',V,T,C) @test_throws DimensionMismatch LAPACK.gemqrt!('R','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('X','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('R','X',V,T,C) C = rand(elty,10,10) V = rand(elty,11,10) T = zeros(elty,10,10) @test_throws DimensionMismatch LAPACK.gemqrt!('R','N',V,T,C) @test_throws DimensionMismatch LAPACK.gemqrt!('L','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('X','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('L','X',V,T,C) # test size(T) = (nb,k) ensures 1 <= nb <= k T = zeros(elty,10,10) V = rand(elty,5,10) @test_throws DimensionMismatch LAPACK.gemqrt!('L','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('X','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('L','X',V,T,C) C = rand(elty,10,10) V = rand(elty,10,10) T = zeros(elty,11,10) @test_throws DimensionMismatch LAPACK.gemqrt!('R','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('X','N',V,T,C) + @test_throws ArgumentError LAPACK.gemqrt!('R','X',V,T,C) @test_throws DimensionMismatch LAPACK.orghr!(1, 10, C, zeros(elty,11)) end @@ -377,8 +437,12 @@ end A = A + transpose(A) #symmetric! B = copy(A) B,ipiv = LAPACK.sytrf!('U',B) + @test_throws ArgumentError LAPACK.sytrf!('X',B) @test triu(inv(A)) ≈ triu(LAPACK.sytri!('U',B,ipiv)) rtol=eps(cond(A)) - @test_throws DimensionMismatch LAPACK.sytrs!('U',B,ipiv,rand(elty,11,5)) + @test_throws ArgumentError LAPACK.sytri!('X',B,ipiv) + temp = rand(elty,11,5) + @test_throws DimensionMismatch LAPACK.sytrs!('U',B,ipiv,temp) + @test_throws ArgumentError LAPACK.sytrs!('X',B,ipiv,temp) @test LAPACK.sytrf!('U',zeros(elty,0,0)) == (zeros(elty,0,0),zeros(BlasInt,0),zero(BlasInt)) end @@ -389,7 +453,10 @@ end B = copy(A) B,ipiv = LAPACK.sytrf_rook!('U', B) @test triu(inv(A)) ≈ triu(LAPACK.sytri_rook!('U', B, ipiv)) rtol=eps(cond(A)) - @test_throws DimensionMismatch LAPACK.sytrs_rook!('U', B, ipiv, rand(elty, 11, 5)) + @test_throws ArgumentError LAPACK.sytri_rook!('X', B, ipiv) + temp = rand(elty, 11, 5) + @test_throws DimensionMismatch LAPACK.sytrs_rook!('U', B, ipiv, temp) + @test_throws ArgumentError LAPACK.sytrs_rook!('X', B, ipiv, temp) @test LAPACK.sytrf_rook!('U',zeros(elty, 0, 0)) == (zeros(elty, 0, 0),zeros(BlasInt, 0),zero(BlasInt)) A = rand(elty, 10, 10) A = A + transpose(A) #symmetric! @@ -398,7 +465,9 @@ end cnd = cond(A) b,A = LAPACK.sysv_rook!('U', A, b) @test b ≈ c rtol=eps(cnd) - @test_throws DimensionMismatch LAPACK.sysv_rook!('U',A,rand(elty,11)) + temp = rand(elty,11) + @test_throws DimensionMismatch LAPACK.sysv_rook!('U',A,temp) + @test_throws ArgumentError LAPACK.sysv_rook!('X',A,temp) # syconvf_rook error handling # way argument is wrong @@ -416,8 +485,11 @@ end A = A + A' #hermitian! B = copy(A) B,ipiv = LAPACK.hetrf!('U',B) - @test_throws DimensionMismatch LAPACK.hetrs!('U',B,ipiv,rand(elty,11,5)) - @test_throws DimensionMismatch LAPACK.hetrs_rook!('U',B,ipiv,rand(elty,11,5)) + temp = rand(elty,11,5) + @test_throws DimensionMismatch LAPACK.hetrs!('U',B,ipiv,temp) + @test_throws ArgumentError LAPACK.hetrs!('X',B,ipiv,temp) + @test_throws DimensionMismatch LAPACK.hetrs_rook!('U',B,ipiv,temp) + @test_throws ArgumentError LAPACK.hetrs_rook!('X',B,ipiv,temp) end end @@ -425,11 +497,21 @@ end @testset for elty in (Float32, Float64) d = rand(elty,10) e = rand(elty,9) - @test_throws DimensionMismatch LAPACK.stev!('U',d,rand(elty,11)) - @test_throws DimensionMismatch LAPACK.stebz!('A','B',zero(elty),zero(elty),0,0,-1.,d,rand(elty,10)) - @test_throws DimensionMismatch LAPACK.stegr!('N','A',d,rand(elty,11),zero(elty),zero(elty),0,0) - @test_throws DimensionMismatch LAPACK.stein!(d,zeros(elty,11),zeros(elty,10),zeros(BlasInt,10),zeros(BlasInt,10)) - @test_throws DimensionMismatch LAPACK.stein!(d,e,zeros(elty,11),zeros(BlasInt,10),zeros(BlasInt,10)) + temp = rand(elty,11) + @test_throws DimensionMismatch LAPACK.stev!('N',d,temp) + @test_throws ArgumentError LAPACK.stev!('X',d,temp) + temp = rand(elty,10) + @test_throws DimensionMismatch LAPACK.stebz!('A','B',zero(elty),zero(elty),0,0,-1.,d,temp) + @test_throws ArgumentError LAPACK.stebz!('X','B',zero(elty),zero(elty),0,0,-1.,d,temp) + @test_throws ArgumentError LAPACK.stebz!('A','X',zero(elty),zero(elty),0,0,-1.,d,temp) + temp11 = rand(elty,11) + @test_throws DimensionMismatch LAPACK.stegr!('N','A',d,temp11,zero(elty),zero(elty),0,0) + @test_throws ArgumentError LAPACK.stegr!('X','A',d,temp11,zero(elty),zero(elty),0,0) + @test_throws ArgumentError LAPACK.stegr!('N','X',d,temp11,zero(elty),zero(elty),0,0) + tempblasint10 = zeros(BlasInt,10) + tempblasint10_2 = zeros(BlasInt,10) + @test_throws DimensionMismatch LAPACK.stein!(d,temp11,temp,tempblasint10,tempblasint10_2) + @test_throws DimensionMismatch LAPACK.stein!(d,e,temp11,tempblasint10,tempblasint10_2) end end @@ -439,7 +521,13 @@ end A = triu(A) B = copy(A) @test inv(A) ≈ LAPACK.trtri!('U','N',B) - @test_throws DimensionMismatch LAPACK.trtrs!('U','N','N',B,zeros(elty,11,10)) + @test_throws ArgumentError LAPACK.trtri!('X','N',B) + @test_throws ArgumentError LAPACK.trtri!('U','X',B) + temp = zeros(elty,11,10) + @test_throws DimensionMismatch LAPACK.trtrs!('U','N','N',B,temp) + @test_throws ArgumentError LAPACK.trtrs!('X','N','N',B,temp) + @test_throws ArgumentError LAPACK.trtrs!('U','X','N',B,temp) + @test_throws ArgumentError LAPACK.trtrs!('U','N','X',B,temp) end end @@ -479,6 +567,8 @@ end LinearAlgebra.LAPACK.larf!('L', v, τ, C1) LinearAlgebra.LAPACK.larf!('R', v, conj(τ), C2) @test C ≈ C2*C1 + + @test_throws ArgumentError LAPACK.larf!('X', v, τ, C1) end end @@ -489,6 +579,8 @@ end @test_throws DimensionMismatch LAPACK.tgsen!(zeros(BlasInt,10),Z,Z,zeros(elty,11,11),Z) @test_throws DimensionMismatch LAPACK.tgsen!(zeros(BlasInt,10),Z,Z,Z,zeros(elty,11,11)) @test_throws DimensionMismatch LAPACK.trsyl!('N','N',Z,Z,zeros(elty,11,11)) + @test_throws ArgumentError LAPACK.trsyl!('X','N',Z,Z,zeros(elty,11,11)) + @test_throws ArgumentError LAPACK.trsyl!('N','X',Z,Z,zeros(elty,11,11)) @test_throws DimensionMismatch LAPACK.tzrzf!(zeros(elty,10,5)) A = triu(rand(elty,4,4)) @@ -508,6 +600,7 @@ end b,A = LAPACK.sysv!('U',A,b) @test b ≈ c @test_throws DimensionMismatch LAPACK.sysv!('U',A,rand(elty,11)) + @test_throws ArgumentError LAPACK.sysv!('X',A,rand(elty,11)) end end @@ -520,14 +613,17 @@ end c = A \ b b,A = LAPACK.hesv!('U',A,b) @test b ≈ c - @test_throws DimensionMismatch LAPACK.hesv!('U',A,rand(elty,11)) + temp = rand(elty,11) + @test_throws DimensionMismatch LAPACK.hesv!('U',A,temp) + @test_throws ArgumentError LAPACK.hesv!('X',A,temp) A = rand(elty,10,10) A = A + A' #hermitian! b = rand(elty,10) c = A \ b b,A = LAPACK.hesv_rook!('U',A,b) @test b ≈ c - @test_throws DimensionMismatch LAPACK.hesv_rook!('U',A,rand(elty,11)) + @test_throws DimensionMismatch LAPACK.hesv_rook!('U',A,temp) + @test_throws ArgumentError LAPACK.hesv_rook!('X',A,temp) end end @@ -563,8 +659,12 @@ end C = copy(B) if elty <: Complex @test A\B ≈ LAPACK.pttrs!('U',rdv,ev,C) - @test_throws DimensionMismatch LAPACK.pttrs!('U',rdv,Vector{elty}(undef,10),C) - @test_throws DimensionMismatch LAPACK.pttrs!('U',rdv,ev,Matrix{elty}(undef,11,11)) + tempvec = Vector{elty}(undef,10) + tempmat = Matrix{elty}(undef,11,11) + @test_throws DimensionMismatch LAPACK.pttrs!('U',rdv,tempvec,C) + @test_throws DimensionMismatch LAPACK.pttrs!('U',rdv,ev,tempmat) + @test_throws ArgumentError LAPACK.pttrs!('X',rdv,tempvec,C) + @test_throws ArgumentError LAPACK.pttrs!('X',rdv,ev,tempmat) else @test A\B ≈ LAPACK.pttrs!(rdv,ev,C) @test_throws DimensionMismatch LAPACK.pttrs!(rdv,Vector{elty}(undef,10),C) @@ -591,6 +691,8 @@ end offsizemat = Matrix{elty}(undef, n+1, n+1) @test_throws DimensionMismatch LAPACK.posv!('U', D, offsizemat) @test_throws DimensionMismatch LAPACK.potrs!('U', D, offsizemat) + @test_throws ArgumentError LAPACK.posv!('X', D, offsizemat) + @test_throws ArgumentError LAPACK.potrs!('X', D, offsizemat) @test LAPACK.potrs!('U',Matrix{elty}(undef,0,0),elty[]) == elty[] end @@ -613,6 +715,10 @@ end B = rand(elty,11,11) @test_throws DimensionMismatch LAPACK.gges!('V','V',A,B) @test_throws DimensionMismatch LAPACK.gges3!('V','V',A,B) + @test_throws ArgumentError LAPACK.gges!('X','V',A,B) + @test_throws ArgumentError LAPACK.gges3!('X','V',A,B) + @test_throws ArgumentError LAPACK.gges!('V','X',A,B) + @test_throws ArgumentError LAPACK.gges3!('V','X',A,B) end end @@ -632,8 +738,14 @@ end select,Vln,Vrn = LAPACK.trevc!('B','S',select,copy(T)) @test Vrn ≈ v @test Vln ≈ Vl - @test_throws ArgumentError LAPACK.trevc!('V','S',select,copy(T)) - @test_throws DimensionMismatch LAPACK.trrfs!('U','N','N',T,rand(elty,10,10),rand(elty,10,11)) + @test_throws ArgumentError LAPACK.trevc!('V','S',select,T) + @test_throws ArgumentError LAPACK.trevc!('R','X',select,T) + temp1010 = rand(elty,10,10) + temp1011 = rand(elty,10,11) + @test_throws DimensionMismatch LAPACK.trrfs!('U','N','N',T,temp1010,temp1011) + @test_throws ArgumentError LAPACK.trrfs!('X','N','N',T,temp1010,temp1011) + @test_throws ArgumentError LAPACK.trrfs!('U','X','N',T,temp1010,temp1011) + @test_throws ArgumentError LAPACK.trrfs!('U','N','X',T,temp1010,temp1011) end end