Annotation of rpl/lapack/lapack/dtrtri.f, revision 1.17

1.8       bertrand    1: *> \brief \b DTRTRI
                      2: *
                      3: *  =========== DOCUMENTATION ===========
                      4: *
1.14      bertrand    5: * Online html documentation available at
                      6: *            http://www.netlib.org/lapack/explore-html/
1.8       bertrand    7: *
                      8: *> \htmlonly
1.14      bertrand    9: *> Download DTRTRI + dependencies
                     10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtrtri.f">
                     11: *> [TGZ]</a>
                     12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtrtri.f">
                     13: *> [ZIP]</a>
                     14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtrtri.f">
1.8       bertrand   15: *> [TXT]</a>
1.14      bertrand   16: *> \endhtmlonly
1.8       bertrand   17: *
                     18: *  Definition:
                     19: *  ===========
                     20: *
                     21: *       SUBROUTINE DTRTRI( UPLO, DIAG, N, A, LDA, INFO )
1.14      bertrand   22: *
1.8       bertrand   23: *       .. Scalar Arguments ..
                     24: *       CHARACTER          DIAG, UPLO
                     25: *       INTEGER            INFO, LDA, N
                     26: *       ..
                     27: *       .. Array Arguments ..
                     28: *       DOUBLE PRECISION   A( LDA, * )
                     29: *       ..
1.14      bertrand   30: *
1.8       bertrand   31: *
                     32: *> \par Purpose:
                     33: *  =============
                     34: *>
                     35: *> \verbatim
                     36: *>
                     37: *> DTRTRI computes the inverse of a real upper or lower triangular
                     38: *> matrix A.
                     39: *>
                     40: *> This is the Level 3 BLAS version of the algorithm.
                     41: *> \endverbatim
                     42: *
                     43: *  Arguments:
                     44: *  ==========
                     45: *
                     46: *> \param[in] UPLO
                     47: *> \verbatim
                     48: *>          UPLO is CHARACTER*1
                     49: *>          = 'U':  A is upper triangular;
                     50: *>          = 'L':  A is lower triangular.
                     51: *> \endverbatim
                     52: *>
                     53: *> \param[in] DIAG
                     54: *> \verbatim
                     55: *>          DIAG is CHARACTER*1
                     56: *>          = 'N':  A is non-unit triangular;
                     57: *>          = 'U':  A is unit triangular.
                     58: *> \endverbatim
                     59: *>
                     60: *> \param[in] N
                     61: *> \verbatim
                     62: *>          N is INTEGER
                     63: *>          The order of the matrix A.  N >= 0.
                     64: *> \endverbatim
                     65: *>
                     66: *> \param[in,out] A
                     67: *> \verbatim
                     68: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
                     69: *>          On entry, the triangular matrix A.  If UPLO = 'U', the
                     70: *>          leading N-by-N upper triangular part of the array A contains
                     71: *>          the upper triangular matrix, and the strictly lower
                     72: *>          triangular part of A is not referenced.  If UPLO = 'L', the
                     73: *>          leading N-by-N lower triangular part of the array A contains
                     74: *>          the lower triangular matrix, and the strictly upper
                     75: *>          triangular part of A is not referenced.  If DIAG = 'U', the
                     76: *>          diagonal elements of A are also not referenced and are
                     77: *>          assumed to be 1.
                     78: *>          On exit, the (triangular) inverse of the original matrix, in
                     79: *>          the same storage format.
                     80: *> \endverbatim
                     81: *>
                     82: *> \param[in] LDA
                     83: *> \verbatim
                     84: *>          LDA is INTEGER
                     85: *>          The leading dimension of the array A.  LDA >= max(1,N).
                     86: *> \endverbatim
                     87: *>
                     88: *> \param[out] INFO
                     89: *> \verbatim
                     90: *>          INFO is INTEGER
                     91: *>          = 0: successful exit
                     92: *>          < 0: if INFO = -i, the i-th argument had an illegal value
                     93: *>          > 0: if INFO = i, A(i,i) is exactly zero.  The triangular
                     94: *>               matrix is singular and its inverse can not be computed.
                     95: *> \endverbatim
                     96: *
                     97: *  Authors:
                     98: *  ========
                     99: *
1.14      bertrand  100: *> \author Univ. of Tennessee
                    101: *> \author Univ. of California Berkeley
                    102: *> \author Univ. of Colorado Denver
                    103: *> \author NAG Ltd.
1.8       bertrand  104: *
                    105: *> \ingroup doubleOTHERcomputational
                    106: *
                    107: *  =====================================================================
1.1       bertrand  108:       SUBROUTINE DTRTRI( UPLO, DIAG, N, A, LDA, INFO )
                    109: *
1.17    ! bertrand  110: *  -- LAPACK computational routine --
1.1       bertrand  111: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    112: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
                    113: *
                    114: *     .. Scalar Arguments ..
                    115:       CHARACTER          DIAG, UPLO
                    116:       INTEGER            INFO, LDA, N
                    117: *     ..
                    118: *     .. Array Arguments ..
                    119:       DOUBLE PRECISION   A( LDA, * )
                    120: *     ..
                    121: *
                    122: *  =====================================================================
                    123: *
                    124: *     .. Parameters ..
                    125:       DOUBLE PRECISION   ONE, ZERO
                    126:       PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
                    127: *     ..
                    128: *     .. Local Scalars ..
                    129:       LOGICAL            NOUNIT, UPPER
                    130:       INTEGER            J, JB, NB, NN
                    131: *     ..
                    132: *     .. External Functions ..
                    133:       LOGICAL            LSAME
                    134:       INTEGER            ILAENV
                    135:       EXTERNAL           LSAME, ILAENV
                    136: *     ..
                    137: *     .. External Subroutines ..
                    138:       EXTERNAL           DTRMM, DTRSM, DTRTI2, XERBLA
                    139: *     ..
                    140: *     .. Intrinsic Functions ..
                    141:       INTRINSIC          MAX, MIN
                    142: *     ..
                    143: *     .. Executable Statements ..
                    144: *
                    145: *     Test the input parameters.
                    146: *
                    147:       INFO = 0
                    148:       UPPER = LSAME( UPLO, 'U' )
                    149:       NOUNIT = LSAME( DIAG, 'N' )
                    150:       IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
                    151:          INFO = -1
                    152:       ELSE IF( .NOT.NOUNIT .AND. .NOT.LSAME( DIAG, 'U' ) ) THEN
                    153:          INFO = -2
                    154:       ELSE IF( N.LT.0 ) THEN
                    155:          INFO = -3
                    156:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
                    157:          INFO = -5
                    158:       END IF
                    159:       IF( INFO.NE.0 ) THEN
                    160:          CALL XERBLA( 'DTRTRI', -INFO )
                    161:          RETURN
                    162:       END IF
                    163: *
                    164: *     Quick return if possible
                    165: *
                    166:       IF( N.EQ.0 )
                    167:      $   RETURN
                    168: *
                    169: *     Check for singularity if non-unit.
                    170: *
                    171:       IF( NOUNIT ) THEN
                    172:          DO 10 INFO = 1, N
                    173:             IF( A( INFO, INFO ).EQ.ZERO )
                    174:      $         RETURN
                    175:    10    CONTINUE
                    176:          INFO = 0
                    177:       END IF
                    178: *
                    179: *     Determine the block size for this environment.
                    180: *
                    181:       NB = ILAENV( 1, 'DTRTRI', UPLO // DIAG, N, -1, -1, -1 )
                    182:       IF( NB.LE.1 .OR. NB.GE.N ) THEN
                    183: *
                    184: *        Use unblocked code
                    185: *
                    186:          CALL DTRTI2( UPLO, DIAG, N, A, LDA, INFO )
                    187:       ELSE
                    188: *
                    189: *        Use blocked code
                    190: *
                    191:          IF( UPPER ) THEN
                    192: *
                    193: *           Compute inverse of upper triangular matrix
                    194: *
                    195:             DO 20 J = 1, N, NB
                    196:                JB = MIN( NB, N-J+1 )
                    197: *
                    198: *              Compute rows 1:j-1 of current block column
                    199: *
                    200:                CALL DTRMM( 'Left', 'Upper', 'No transpose', DIAG, J-1,
                    201:      $                     JB, ONE, A, LDA, A( 1, J ), LDA )
                    202:                CALL DTRSM( 'Right', 'Upper', 'No transpose', DIAG, J-1,
                    203:      $                     JB, -ONE, A( J, J ), LDA, A( 1, J ), LDA )
                    204: *
                    205: *              Compute inverse of current diagonal block
                    206: *
                    207:                CALL DTRTI2( 'Upper', DIAG, JB, A( J, J ), LDA, INFO )
                    208:    20       CONTINUE
                    209:          ELSE
                    210: *
                    211: *           Compute inverse of lower triangular matrix
                    212: *
                    213:             NN = ( ( N-1 ) / NB )*NB + 1
                    214:             DO 30 J = NN, 1, -NB
                    215:                JB = MIN( NB, N-J+1 )
                    216:                IF( J+JB.LE.N ) THEN
                    217: *
                    218: *                 Compute rows j+jb:n of current block column
                    219: *
                    220:                   CALL DTRMM( 'Left', 'Lower', 'No transpose', DIAG,
                    221:      $                        N-J-JB+1, JB, ONE, A( J+JB, J+JB ), LDA,
                    222:      $                        A( J+JB, J ), LDA )
                    223:                   CALL DTRSM( 'Right', 'Lower', 'No transpose', DIAG,
                    224:      $                        N-J-JB+1, JB, -ONE, A( J, J ), LDA,
                    225:      $                        A( J+JB, J ), LDA )
                    226:                END IF
                    227: *
                    228: *              Compute inverse of current diagonal block
                    229: *
                    230:                CALL DTRTI2( 'Lower', DIAG, JB, A( J, J ), LDA, INFO )
                    231:    30       CONTINUE
                    232:          END IF
                    233:       END IF
                    234: *
                    235:       RETURN
                    236: *
                    237: *     End of DTRTRI
                    238: *
                    239:       END

CVSweb interface <joel.bertrand@systella.fr>