Annotation of rpl/lapack/lapack/dgetri.f, revision 1.1

1.1     ! bertrand    1:       SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO )
        !             2: *
        !             3: *  -- LAPACK routine (version 3.2) --
        !             4: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
        !             5: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
        !             6: *     November 2006
        !             7: *
        !             8: *     .. Scalar Arguments ..
        !             9:       INTEGER            INFO, LDA, LWORK, N
        !            10: *     ..
        !            11: *     .. Array Arguments ..
        !            12:       INTEGER            IPIV( * )
        !            13:       DOUBLE PRECISION   A( LDA, * ), WORK( * )
        !            14: *     ..
        !            15: *
        !            16: *  Purpose
        !            17: *  =======
        !            18: *
        !            19: *  DGETRI computes the inverse of a matrix using the LU factorization
        !            20: *  computed by DGETRF.
        !            21: *
        !            22: *  This method inverts U and then computes inv(A) by solving the system
        !            23: *  inv(A)*L = inv(U) for inv(A).
        !            24: *
        !            25: *  Arguments
        !            26: *  =========
        !            27: *
        !            28: *  N       (input) INTEGER
        !            29: *          The order of the matrix A.  N >= 0.
        !            30: *
        !            31: *  A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
        !            32: *          On entry, the factors L and U from the factorization
        !            33: *          A = P*L*U as computed by DGETRF.
        !            34: *          On exit, if INFO = 0, the inverse of the original matrix A.
        !            35: *
        !            36: *  LDA     (input) INTEGER
        !            37: *          The leading dimension of the array A.  LDA >= max(1,N).
        !            38: *
        !            39: *  IPIV    (input) INTEGER array, dimension (N)
        !            40: *          The pivot indices from DGETRF; for 1<=i<=N, row i of the
        !            41: *          matrix was interchanged with row IPIV(i).
        !            42: *
        !            43: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
        !            44: *          On exit, if INFO=0, then WORK(1) returns the optimal LWORK.
        !            45: *
        !            46: *  LWORK   (input) INTEGER
        !            47: *          The dimension of the array WORK.  LWORK >= max(1,N).
        !            48: *          For optimal performance LWORK >= N*NB, where NB is
        !            49: *          the optimal blocksize returned by ILAENV.
        !            50: *
        !            51: *          If LWORK = -1, then a workspace query is assumed; the routine
        !            52: *          only calculates the optimal size of the WORK array, returns
        !            53: *          this value as the first entry of the WORK array, and no error
        !            54: *          message related to LWORK is issued by XERBLA.
        !            55: *
        !            56: *  INFO    (output) INTEGER
        !            57: *          = 0:  successful exit
        !            58: *          < 0:  if INFO = -i, the i-th argument had an illegal value
        !            59: *          > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
        !            60: *                singular and its inverse could not be computed.
        !            61: *
        !            62: *  =====================================================================
        !            63: *
        !            64: *     .. Parameters ..
        !            65:       DOUBLE PRECISION   ZERO, ONE
        !            66:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
        !            67: *     ..
        !            68: *     .. Local Scalars ..
        !            69:       LOGICAL            LQUERY
        !            70:       INTEGER            I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB,
        !            71:      $                   NBMIN, NN
        !            72: *     ..
        !            73: *     .. External Functions ..
        !            74:       INTEGER            ILAENV
        !            75:       EXTERNAL           ILAENV
        !            76: *     ..
        !            77: *     .. External Subroutines ..
        !            78:       EXTERNAL           DGEMM, DGEMV, DSWAP, DTRSM, DTRTRI, XERBLA
        !            79: *     ..
        !            80: *     .. Intrinsic Functions ..
        !            81:       INTRINSIC          MAX, MIN
        !            82: *     ..
        !            83: *     .. Executable Statements ..
        !            84: *
        !            85: *     Test the input parameters.
        !            86: *
        !            87:       INFO = 0
        !            88:       NB = ILAENV( 1, 'DGETRI', ' ', N, -1, -1, -1 )
        !            89:       LWKOPT = N*NB
        !            90:       WORK( 1 ) = LWKOPT
        !            91:       LQUERY = ( LWORK.EQ.-1 )
        !            92:       IF( N.LT.0 ) THEN
        !            93:          INFO = -1
        !            94:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
        !            95:          INFO = -3
        !            96:       ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN
        !            97:          INFO = -6
        !            98:       END IF
        !            99:       IF( INFO.NE.0 ) THEN
        !           100:          CALL XERBLA( 'DGETRI', -INFO )
        !           101:          RETURN
        !           102:       ELSE IF( LQUERY ) THEN
        !           103:          RETURN
        !           104:       END IF
        !           105: *
        !           106: *     Quick return if possible
        !           107: *
        !           108:       IF( N.EQ.0 )
        !           109:      $   RETURN
        !           110: *
        !           111: *     Form inv(U).  If INFO > 0 from DTRTRI, then U is singular,
        !           112: *     and the inverse is not computed.
        !           113: *
        !           114:       CALL DTRTRI( 'Upper', 'Non-unit', N, A, LDA, INFO )
        !           115:       IF( INFO.GT.0 )
        !           116:      $   RETURN
        !           117: *
        !           118:       NBMIN = 2
        !           119:       LDWORK = N
        !           120:       IF( NB.GT.1 .AND. NB.LT.N ) THEN
        !           121:          IWS = MAX( LDWORK*NB, 1 )
        !           122:          IF( LWORK.LT.IWS ) THEN
        !           123:             NB = LWORK / LDWORK
        !           124:             NBMIN = MAX( 2, ILAENV( 2, 'DGETRI', ' ', N, -1, -1, -1 ) )
        !           125:          END IF
        !           126:       ELSE
        !           127:          IWS = N
        !           128:       END IF
        !           129: *
        !           130: *     Solve the equation inv(A)*L = inv(U) for inv(A).
        !           131: *
        !           132:       IF( NB.LT.NBMIN .OR. NB.GE.N ) THEN
        !           133: *
        !           134: *        Use unblocked code.
        !           135: *
        !           136:          DO 20 J = N, 1, -1
        !           137: *
        !           138: *           Copy current column of L to WORK and replace with zeros.
        !           139: *
        !           140:             DO 10 I = J + 1, N
        !           141:                WORK( I ) = A( I, J )
        !           142:                A( I, J ) = ZERO
        !           143:    10       CONTINUE
        !           144: *
        !           145: *           Compute current column of inv(A).
        !           146: *
        !           147:             IF( J.LT.N )
        !           148:      $         CALL DGEMV( 'No transpose', N, N-J, -ONE, A( 1, J+1 ),
        !           149:      $                     LDA, WORK( J+1 ), 1, ONE, A( 1, J ), 1 )
        !           150:    20    CONTINUE
        !           151:       ELSE
        !           152: *
        !           153: *        Use blocked code.
        !           154: *
        !           155:          NN = ( ( N-1 ) / NB )*NB + 1
        !           156:          DO 50 J = NN, 1, -NB
        !           157:             JB = MIN( NB, N-J+1 )
        !           158: *
        !           159: *           Copy current block column of L to WORK and replace with
        !           160: *           zeros.
        !           161: *
        !           162:             DO 40 JJ = J, J + JB - 1
        !           163:                DO 30 I = JJ + 1, N
        !           164:                   WORK( I+( JJ-J )*LDWORK ) = A( I, JJ )
        !           165:                   A( I, JJ ) = ZERO
        !           166:    30          CONTINUE
        !           167:    40       CONTINUE
        !           168: *
        !           169: *           Compute current block column of inv(A).
        !           170: *
        !           171:             IF( J+JB.LE.N )
        !           172:      $         CALL DGEMM( 'No transpose', 'No transpose', N, JB,
        !           173:      $                     N-J-JB+1, -ONE, A( 1, J+JB ), LDA,
        !           174:      $                     WORK( J+JB ), LDWORK, ONE, A( 1, J ), LDA )
        !           175:             CALL DTRSM( 'Right', 'Lower', 'No transpose', 'Unit', N, JB,
        !           176:      $                  ONE, WORK( J ), LDWORK, A( 1, J ), LDA )
        !           177:    50    CONTINUE
        !           178:       END IF
        !           179: *
        !           180: *     Apply column interchanges.
        !           181: *
        !           182:       DO 60 J = N - 1, 1, -1
        !           183:          JP = IPIV( J )
        !           184:          IF( JP.NE.J )
        !           185:      $      CALL DSWAP( N, A( 1, J ), 1, A( 1, JP ), 1 )
        !           186:    60 CONTINUE
        !           187: *
        !           188:       WORK( 1 ) = IWS
        !           189:       RETURN
        !           190: *
        !           191: *     End of DGETRI
        !           192: *
        !           193:       END

CVSweb interface <joel.bertrand@systella.fr>