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

NAME

       la_gbrcond - la_gbrcond: Skeel condition number estimate

SYNOPSIS

   Functions
       real function cla_gbrcond_c (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, c, capply,
           info, work, rwork)
           CLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for
           general banded matrices.
       real function cla_gbrcond_x (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, x, info, work,
           rwork)
           CLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general
           banded matrices.
       double precision function dla_gbrcond (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv,
           cmode, c, info, work, iwork)
           DLA_GBRCOND estimates the Skeel condition number for a general banded matrix.
       real function sla_gbrcond (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, cmode, c, info,
           work, iwork)
           SLA_GBRCOND estimates the Skeel condition number for a general banded matrix.
       double precision function zla_gbrcond_c (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, c,
           capply, info, work, rwork)
           ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for
           general banded matrices.
       double precision function zla_gbrcond_x (trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, x,
           info, work, rwork)
           ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general
           banded matrices.

Detailed Description

Function Documentation

   real function cla_gbrcond_c (character trans, integer n, integer kl, integer ku, complex,
       dimension( ldab, * ) ab, integer ldab, complex, dimension( ldafb, * ) afb, integer ldafb,
       integer, dimension( * ) ipiv, real, dimension( * ) c, logical capply, integer info,
       complex, dimension( * ) work, real, dimension( * ) rwork)
       CLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for
       general banded matrices.

       Purpose:

               CLA_GBRCOND_C Computes the infinity norm condition number of
               op(A) * inv(diag(C)) where C is a REAL vector.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is COMPLEX array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is COMPLEX array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by CGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by CGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           C

                     C is REAL array, dimension (N)
                The vector C in the formula op(A) * inv(diag(C)).

           CAPPLY

                     CAPPLY is LOGICAL
                If .TRUE. then access the vector C in the formula above.

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

                     WORK is COMPLEX array, dimension (2*N).
                Workspace.

           RWORK

                     RWORK is REAL array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   real function cla_gbrcond_x (character trans, integer n, integer kl, integer ku, complex,
       dimension( ldab, * ) ab, integer ldab, complex, dimension( ldafb, * ) afb, integer ldafb,
       integer, dimension( * ) ipiv, complex, dimension( * ) x, integer info, complex, dimension(
       * ) work, real, dimension( * ) rwork)
       CLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general
       banded matrices.

       Purpose:

               CLA_GBRCOND_X Computes the infinity norm condition number of
               op(A) * diag(X) where X is a COMPLEX vector.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is COMPLEX array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is COMPLEX array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by CGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by CGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           X

                     X is COMPLEX array, dimension (N)
                The vector X in the formula op(A) * diag(X).

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

                     WORK is COMPLEX array, dimension (2*N).
                Workspace.

           RWORK

                     RWORK is REAL array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   double precision function dla_gbrcond (character trans, integer n, integer kl, integer ku,
       double precision, dimension( ldab, * ) ab, integer ldab, double precision, dimension(
       ldafb, * ) afb, integer ldafb, integer, dimension( * ) ipiv, integer cmode, double
       precision, dimension( * ) c, integer info, double precision, dimension( * ) work, integer,
       dimension( * ) iwork)
       DLA_GBRCOND estimates the Skeel condition number for a general banded matrix.

       Purpose:

               DLA_GBRCOND Estimates the Skeel condition number of  op(A) * op2(C)
               where op2 is determined by CMODE as follows
               CMODE =  1    op2(C) = C
               CMODE =  0    op2(C) = I
               CMODE = -1    op2(C) = inv(C)
               The Skeel condition number  cond(A) = norminf( |inv(A)||A| )
               is computed by computing scaling factors R such that
               diag(R)*A*op2(C) is row equilibrated and computing the standard
               infinity-norm condition number.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is DOUBLE PRECISION array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is DOUBLE PRECISION array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by DGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by DGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           CMODE

                     CMODE is INTEGER
                Determines op2(C) in the formula op(A) * op2(C) as follows:
                CMODE =  1    op2(C) = C
                CMODE =  0    op2(C) = I
                CMODE = -1    op2(C) = inv(C)

           C

                     C is DOUBLE PRECISION array, dimension (N)
                The vector C in the formula op(A) * op2(C).

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

                     WORK is DOUBLE PRECISION array, dimension (5*N).
                Workspace.

           IWORK

                     IWORK is INTEGER array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   real function sla_gbrcond (character trans, integer n, integer kl, integer ku, real,
       dimension( ldab, * ) ab, integer ldab, real, dimension( ldafb, * ) afb, integer ldafb,
       integer, dimension( * ) ipiv, integer cmode, real, dimension( * ) c, integer info, real,
       dimension( * ) work, integer, dimension( * ) iwork)
       SLA_GBRCOND estimates the Skeel condition number for a general banded matrix.

       Purpose:

               SLA_GBRCOND Estimates the Skeel condition number of  op(A) * op2(C)
               where op2 is determined by CMODE as follows
               CMODE =  1    op2(C) = C
               CMODE =  0    op2(C) = I
               CMODE = -1    op2(C) = inv(C)
               The Skeel condition number  cond(A) = norminf( |inv(A)||A| )
               is computed by computing scaling factors R such that
               diag(R)*A*op2(C) is row equilibrated and computing the standard
               infinity-norm condition number.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is REAL array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is REAL array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by SGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by SGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           CMODE

                     CMODE is INTEGER
                Determines op2(C) in the formula op(A) * op2(C) as follows:
                CMODE =  1    op2(C) = C
                CMODE =  0    op2(C) = I
                CMODE = -1    op2(C) = inv(C)

           C

                     C is REAL array, dimension (N)
                The vector C in the formula op(A) * op2(C).

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

                     WORK is REAL array, dimension (5*N).
                Workspace.

           IWORK

                     IWORK is INTEGER array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   double precision function zla_gbrcond_c (character trans, integer n, integer kl, integer ku,
       complex*16, dimension( ldab, * ) ab, integer ldab, complex*16, dimension( ldafb, * ) afb,
       integer ldafb, integer, dimension( * ) ipiv, double precision, dimension( * ) c, logical
       capply, integer info, complex*16, dimension( * ) work, double precision, dimension( * )
       rwork)
       ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for
       general banded matrices.

       Purpose:

               ZLA_GBRCOND_C Computes the infinity norm condition number of
               op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is COMPLEX*16 array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is COMPLEX*16 array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by ZGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by ZGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           C

                     C is DOUBLE PRECISION array, dimension (N)
                The vector C in the formula op(A) * inv(diag(C)).

           CAPPLY

                     CAPPLY is LOGICAL
                If .TRUE. then access the vector C in the formula above.

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

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

           RWORK

                     RWORK is DOUBLE PRECISION array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   double precision function zla_gbrcond_x (character trans, integer n, integer kl, integer ku,
       complex*16, dimension( ldab, * ) ab, integer ldab, complex*16, dimension( ldafb, * ) afb,
       integer ldafb, integer, dimension( * ) ipiv, complex*16, dimension( * ) x, integer info,
       complex*16, dimension( * ) work, double precision, dimension( * ) rwork)
       ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general
       banded matrices.

       Purpose:

               ZLA_GBRCOND_X Computes the infinity norm condition number of
               op(A) * diag(X) where X is a COMPLEX*16 vector.

       Parameters
           TRANS

                     TRANS is CHARACTER*1
                Specifies the form of the system of equations:
                  = 'N':  A * X = B     (No transpose)
                  = 'T':  A**T * X = B  (Transpose)
                  = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)

           N

                     N is INTEGER
                The number of linear equations, i.e., the order of the
                matrix A.  N >= 0.

           KL

                     KL is INTEGER
                The number of subdiagonals within the band of A.  KL >= 0.

           KU

                     KU is INTEGER
                The number of superdiagonals within the band of A.  KU >= 0.

           AB

                     AB is COMPLEX*16 array, dimension (LDAB,N)
                On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                The j-th column of A is stored in the j-th column of the
                array AB as follows:
                AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

           LDAB

                     LDAB is INTEGER
                The leading dimension of the array AB.  LDAB >= KL+KU+1.

           AFB

                     AFB is COMPLEX*16 array, dimension (LDAFB,N)
                Details of the LU factorization of the band matrix A, as
                computed by ZGBTRF.  U is stored as an upper triangular
                band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                and the multipliers used during the factorization are stored
                in rows KL+KU+2 to 2*KL+KU+1.

           LDAFB

                     LDAFB is INTEGER
                The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.

           IPIV

                     IPIV is INTEGER array, dimension (N)
                The pivot indices from the factorization A = P*L*U
                as computed by ZGBTRF; row i of the matrix was interchanged
                with row IPIV(i).

           X

                     X is COMPLEX*16 array, dimension (N)
                The vector X in the formula op(A) * diag(X).

           INFO

                     INFO is INTEGER
                  = 0:  Successful exit.
                i > 0:  The ith argument is invalid.

           WORK

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

           RWORK

                     RWORK is DOUBLE PRECISION array, dimension (N).
                Workspace.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

Author

       Generated automatically by Doxygen for LAPACK from the source code.