Provided by: liblapack-doc_3.12.0-3build1.1_all bug

NAME

       getri - getri: triangular inverse

SYNOPSIS

   Functions
       subroutine cgetri (n, a, lda, ipiv, work, lwork, info)
           CGETRI
       subroutine dgetri (n, a, lda, ipiv, work, lwork, info)
           DGETRI
       subroutine sgetri (n, a, lda, ipiv, work, lwork, info)
           SGETRI
       subroutine zgetri (n, a, lda, ipiv, work, lwork, info)
           ZGETRI

Detailed Description

Function Documentation

   subroutine cgetri (integer n, complex, dimension( lda, * ) a, integer lda, integer, dimension(
       * ) ipiv, complex, dimension( * ) work, integer lwork, integer info)
       CGETRI

       Purpose:

            CGETRI computes the inverse of a matrix using the LU factorization
            computed by CGETRF.

            This method inverts U and then computes inv(A) by solving the system
            inv(A)*L = inv(U) for inv(A).

       Parameters
           N

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

           A

                     A is COMPLEX array, dimension (LDA,N)
                     On entry, the factors L and U from the factorization
                     A = P*L*U as computed by CGETRF.
                     On exit, if INFO = 0, the inverse of the original matrix A.

           LDA

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

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     The pivot indices from CGETRF; for 1<=i<=N, row i of the
                     matrix was interchanged with row IPIV(i).

           WORK

                     WORK is COMPLEX array, dimension (MAX(1,LWORK))
                     On exit, if INFO=0, then WORK(1) returns the optimal LWORK.

           LWORK

                     LWORK is INTEGER
                     The dimension of the array WORK.  LWORK >= max(1,N).
                     For optimal performance LWORK >= N*NB, where NB is
                     the optimal blocksize returned by ILAENV.

                     If LWORK = -1, then a workspace query is assumed; the routine
                     only calculates the optimal size of the WORK array, returns
                     this value as the first entry of the WORK array, and no error
                     message related to LWORK is issued by XERBLA.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
                           singular and its inverse could not be computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   subroutine dgetri (integer n, double precision, dimension( lda, * ) a, integer lda, integer,
       dimension( * ) ipiv, double precision, dimension( * ) work, integer lwork, integer info)
       DGETRI

       Purpose:

            DGETRI computes the inverse of a matrix using the LU factorization
            computed by DGETRF.

            This method inverts U and then computes inv(A) by solving the system
            inv(A)*L = inv(U) for inv(A).

       Parameters
           N

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

           A

                     A is DOUBLE PRECISION array, dimension (LDA,N)
                     On entry, the factors L and U from the factorization
                     A = P*L*U as computed by DGETRF.
                     On exit, if INFO = 0, the inverse of the original matrix A.

           LDA

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

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     The pivot indices from DGETRF; for 1<=i<=N, row i of the
                     matrix was interchanged with row IPIV(i).

           WORK

                     WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
                     On exit, if INFO=0, then WORK(1) returns the optimal LWORK.

           LWORK

                     LWORK is INTEGER
                     The dimension of the array WORK.  LWORK >= max(1,N).
                     For optimal performance LWORK >= N*NB, where NB is
                     the optimal blocksize returned by ILAENV.

                     If LWORK = -1, then a workspace query is assumed; the routine
                     only calculates the optimal size of the WORK array, returns
                     this value as the first entry of the WORK array, and no error
                     message related to LWORK is issued by XERBLA.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
                           singular and its inverse could not be computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   subroutine sgetri (integer n, real, dimension( lda, * ) a, integer lda, integer, dimension( *
       ) ipiv, real, dimension( * ) work, integer lwork, integer info)
       SGETRI

       Purpose:

            SGETRI computes the inverse of a matrix using the LU factorization
            computed by SGETRF.

            This method inverts U and then computes inv(A) by solving the system
            inv(A)*L = inv(U) for inv(A).

       Parameters
           N

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

           A

                     A is REAL array, dimension (LDA,N)
                     On entry, the factors L and U from the factorization
                     A = P*L*U as computed by SGETRF.
                     On exit, if INFO = 0, the inverse of the original matrix A.

           LDA

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

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     The pivot indices from SGETRF; for 1<=i<=N, row i of the
                     matrix was interchanged with row IPIV(i).

           WORK

                     WORK is REAL array, dimension (MAX(1,LWORK))
                     On exit, if INFO=0, then WORK(1) returns the optimal LWORK.

           LWORK

                     LWORK is INTEGER
                     The dimension of the array WORK.  LWORK >= max(1,N).
                     For optimal performance LWORK >= N*NB, where NB is
                     the optimal blocksize returned by ILAENV.

                     If LWORK = -1, then a workspace query is assumed; the routine
                     only calculates the optimal size of the WORK array, returns
                     this value as the first entry of the WORK array, and no error
                     message related to LWORK is issued by XERBLA.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
                           singular and its inverse could not be computed.

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

   subroutine zgetri (integer n, complex*16, dimension( lda, * ) a, integer lda, integer,
       dimension( * ) ipiv, complex*16, dimension( * ) work, integer lwork, integer info)
       ZGETRI

       Purpose:

            ZGETRI computes the inverse of a matrix using the LU factorization
            computed by ZGETRF.

            This method inverts U and then computes inv(A) by solving the system
            inv(A)*L = inv(U) for inv(A).

       Parameters
           N

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

           A

                     A is COMPLEX*16 array, dimension (LDA,N)
                     On entry, the factors L and U from the factorization
                     A = P*L*U as computed by ZGETRF.
                     On exit, if INFO = 0, the inverse of the original matrix A.

           LDA

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

           IPIV

                     IPIV is INTEGER array, dimension (N)
                     The pivot indices from ZGETRF; for 1<=i<=N, row i of the
                     matrix was interchanged with row IPIV(i).

           WORK

                     WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
                     On exit, if INFO=0, then WORK(1) returns the optimal LWORK.

           LWORK

                     LWORK is INTEGER
                     The dimension of the array WORK.  LWORK >= max(1,N).
                     For optimal performance LWORK >= N*NB, where NB is
                     the optimal blocksize returned by ILAENV.

                     If LWORK = -1, then a workspace query is assumed; the routine
                     only calculates the optimal size of the WORK array, returns
                     this value as the first entry of the WORK array, and no error
                     message related to LWORK is issued by XERBLA.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value
                     > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
                           singular and its inverse could not be computed.

       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.