Provided by: liblapack-doc_3.12.0-3build1_all bug

NAME

       hecon_3 - {he,sy}con_3: condition number estimate

SYNOPSIS

   Functions
       subroutine checon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
           CHECON_3
       subroutine csycon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
           CSYCON_3
       subroutine dsycon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, iwork, info)
           DSYCON_3
       subroutine ssycon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, iwork, info)
           SSYCON_3
       subroutine zhecon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
           ZHECON_3
       subroutine zsycon_3 (uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
           ZSYCON_3

Detailed Description

Function Documentation

   subroutine checon_3 (character uplo, integer n, complex, dimension( lda, * ) a, integer lda,
       complex, dimension( * ) e, integer, dimension( * ) ipiv, real anorm, real rcond, complex,
       dimension( * ) work, integer info)
       CHECON_3

       Purpose:

            CHECON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a complex Hermitian matrix A using the factorization
            computed by CHETRF_RK or CHETRF_BK:

               A = P*U*D*(U**H)*(P**T) or A = P*L*D*(L**H)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**H (or L**H) is the conjugate of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is Hermitian and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver CHETRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**H)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**H)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is COMPLEX array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by CHETRF_RK and CHETRF_BK:
                       a) ONLY diagonal elements of the Hermitian block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is COMPLEX array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the Hermitian block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by CHETRF_RK or CHETRF_BK.

           ANORM

                     ANORM is REAL
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is REAL
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is COMPLEX array, dimension (2*N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

   subroutine csycon_3 (character uplo, integer n, complex, dimension( lda, * ) a, integer lda,
       complex, dimension( * ) e, integer, dimension( * ) ipiv, real anorm, real rcond, complex,
       dimension( * ) work, integer info)
       CSYCON_3

       Purpose:

            CSYCON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a complex symmetric matrix A using the factorization
            computed by CSYTRF_RK or CSYTRF_BK:

               A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**T (or L**T) is the transpose of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is symmetric and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver CSYTRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**T)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**T)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is COMPLEX array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by CSYTRF_RK and CSYTRF_BK:
                       a) ONLY diagonal elements of the symmetric block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is COMPLEX array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the symmetric block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by CSYTRF_RK or CSYTRF_BK.

           ANORM

                     ANORM is REAL
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is REAL
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is COMPLEX array, dimension (2*N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

   subroutine dsycon_3 (character uplo, integer n, double precision, dimension( lda, * ) a,
       integer lda, double precision, dimension( * ) e, integer, dimension( * ) ipiv, double
       precision anorm, double precision rcond, double precision, dimension( * ) work, integer,
       dimension( * ) iwork, integer info)
       DSYCON_3

       Purpose:

            DSYCON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a real symmetric matrix A using the factorization
            computed by DSYTRF_RK or DSYTRF_BK:

               A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**T (or L**T) is the transpose of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is symmetric and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver DSYTRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**T)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**T)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is DOUBLE PRECISION array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by DSYTRF_RK and DSYTRF_BK:
                       a) ONLY diagonal elements of the symmetric block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is DOUBLE PRECISION array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the symmetric block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by DSYTRF_RK or DSYTRF_BK.

           ANORM

                     ANORM is DOUBLE PRECISION
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is DOUBLE PRECISION
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is DOUBLE PRECISION array, dimension (2*N)

           IWORK

                     IWORK is INTEGER array, dimension (N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

   subroutine ssycon_3 (character uplo, integer n, real, dimension( lda, * ) a, integer lda,
       real, dimension( * ) e, integer, dimension( * ) ipiv, real anorm, real rcond, real,
       dimension( * ) work, integer, dimension( * ) iwork, integer info)
       SSYCON_3

       Purpose:

            SSYCON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a real symmetric matrix A using the factorization
            computed by DSYTRF_RK or DSYTRF_BK:

               A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**T (or L**T) is the transpose of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is symmetric and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver SSYTRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**T)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**T)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is REAL array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by SSYTRF_RK and SSYTRF_BK:
                       a) ONLY diagonal elements of the symmetric block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is REAL array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the symmetric block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by SSYTRF_RK or SSYTRF_BK.

           ANORM

                     ANORM is REAL
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is REAL
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is REAL array, dimension (2*N)

           IWORK

                     IWORK is INTEGER array, dimension (N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

   subroutine zhecon_3 (character uplo, integer n, complex*16, dimension( lda, * ) a, integer
       lda, complex*16, dimension( * ) e, integer, dimension( * ) ipiv, double precision anorm,
       double precision rcond, complex*16, dimension( * ) work, integer info)
       ZHECON_3

       Purpose:

            ZHECON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a complex Hermitian matrix A using the factorization
            computed by ZHETRF_RK or ZHETRF_BK:

               A = P*U*D*(U**H)*(P**T) or A = P*L*D*(L**H)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**H (or L**H) is the conjugate of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is Hermitian and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver ZHETRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**H)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**H)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is COMPLEX*16 array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by ZHETRF_RK and ZHETRF_BK:
                       a) ONLY diagonal elements of the Hermitian block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is COMPLEX*16 array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the Hermitian block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by ZHETRF_RK or ZHETRF_BK.

           ANORM

                     ANORM is DOUBLE PRECISION
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is DOUBLE PRECISION
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is COMPLEX*16 array, dimension (2*N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

   subroutine zsycon_3 (character uplo, integer n, complex*16, dimension( lda, * ) a, integer
       lda, complex*16, dimension( * ) e, integer, dimension( * ) ipiv, double precision anorm,
       double precision rcond, complex*16, dimension( * ) work, integer info)
       ZSYCON_3

       Purpose:

            ZSYCON_3 estimates the reciprocal of the condition number (in the
            1-norm) of a complex symmetric matrix A using the factorization
            computed by ZSYTRF_RK or ZSYTRF_BK:

               A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),

            where U (or L) is unit upper (or lower) triangular matrix,
            U**T (or L**T) is the transpose of U (or L), P is a permutation
            matrix, P**T is the transpose of P, and D is symmetric and block
            diagonal with 1-by-1 and 2-by-2 diagonal blocks.

            An estimate is obtained for norm(inv(A)), and the reciprocal of the
            condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
            This routine uses BLAS3 solver ZSYTRS_3.

       Parameters
           UPLO

                     UPLO is CHARACTER*1
                     Specifies whether the details of the factorization are
                     stored as an upper or lower triangular matrix:
                     = 'U':  Upper triangular, form is A = P*U*D*(U**T)*(P**T);
                     = 'L':  Lower triangular, form is A = P*L*D*(L**T)*(P**T).

           N

                     N is INTEGER
                     The order of the matrix A.  N >= 0.

           A

                     A is COMPLEX*16 array, dimension (LDA,N)
                     Diagonal of the block diagonal matrix D and factors U or L
                     as computed by ZSYTRF_RK and ZSYTRF_BK:
                       a) ONLY diagonal elements of the symmetric block diagonal
                          matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
                          (superdiagonal (or subdiagonal) elements of D
                           should be provided on entry in array E), and
                       b) If UPLO = 'U': factor U in the superdiagonal part of A.
                          If UPLO = 'L': factor L in the subdiagonal part of A.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,N).

           E

                     E is COMPLEX*16 array, dimension (N)
                     On entry, contains the superdiagonal (or subdiagonal)
                     elements of the symmetric block diagonal matrix D
                     with 1-by-1 or 2-by-2 diagonal blocks, where
                     If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not referenced;
                     If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

                     NOTE: For 1-by-1 diagonal block D(k), where
                     1 <= k <= N, the element E(k) is not referenced in both
                     UPLO = 'U' or UPLO = 'L' cases.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     Details of the interchanges and the block structure of D
                     as determined by ZSYTRF_RK or ZSYTRF_BK.

           ANORM

                     ANORM is DOUBLE PRECISION
                     The 1-norm of the original matrix A.

           RCOND

                     RCOND is DOUBLE PRECISION
                     The reciprocal of the condition number of the matrix A,
                     computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
                     estimate of the 1-norm of inv(A) computed in this routine.

           WORK

                     WORK is COMPLEX*16 array, dimension (2*N)

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

             June 2017,  Igor Kozachenko,
                             Computer Science Division,
                             University of California, Berkeley

             September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                             School of Mathematics,
                             University of Manchester

Author

       Generated automatically by Doxygen for LAPACK from the source code.