LAPACK 3.3.0
|
00001 SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO ) 00002 * 00003 * -- LAPACK routine (version 3.2) -- 00004 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00005 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00006 * November 2006 00007 * 00008 * .. Scalar Arguments .. 00009 INTEGER INFO, LDA, LWORK, N 00010 * .. 00011 * .. Array Arguments .. 00012 INTEGER IPIV( * ) 00013 DOUBLE PRECISION A( LDA, * ), WORK( * ) 00014 * .. 00015 * 00016 * Purpose 00017 * ======= 00018 * 00019 * DGETRI computes the inverse of a matrix using the LU factorization 00020 * computed by DGETRF. 00021 * 00022 * This method inverts U and then computes inv(A) by solving the system 00023 * inv(A)*L = inv(U) for inv(A). 00024 * 00025 * Arguments 00026 * ========= 00027 * 00028 * N (input) INTEGER 00029 * The order of the matrix A. N >= 0. 00030 * 00031 * A (input/output) DOUBLE PRECISION array, dimension (LDA,N) 00032 * On entry, the factors L and U from the factorization 00033 * A = P*L*U as computed by DGETRF. 00034 * On exit, if INFO = 0, the inverse of the original matrix A. 00035 * 00036 * LDA (input) INTEGER 00037 * The leading dimension of the array A. LDA >= max(1,N). 00038 * 00039 * IPIV (input) INTEGER array, dimension (N) 00040 * The pivot indices from DGETRF; for 1<=i<=N, row i of the 00041 * matrix was interchanged with row IPIV(i). 00042 * 00043 * WORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) 00044 * On exit, if INFO=0, then WORK(1) returns the optimal LWORK. 00045 * 00046 * LWORK (input) INTEGER 00047 * The dimension of the array WORK. LWORK >= max(1,N). 00048 * For optimal performance LWORK >= N*NB, where NB is 00049 * the optimal blocksize returned by ILAENV. 00050 * 00051 * If LWORK = -1, then a workspace query is assumed; the routine 00052 * only calculates the optimal size of the WORK array, returns 00053 * this value as the first entry of the WORK array, and no error 00054 * message related to LWORK is issued by XERBLA. 00055 * 00056 * INFO (output) INTEGER 00057 * = 0: successful exit 00058 * < 0: if INFO = -i, the i-th argument had an illegal value 00059 * > 0: if INFO = i, U(i,i) is exactly zero; the matrix is 00060 * singular and its inverse could not be computed. 00061 * 00062 * ===================================================================== 00063 * 00064 * .. Parameters .. 00065 DOUBLE PRECISION ZERO, ONE 00066 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00067 * .. 00068 * .. Local Scalars .. 00069 LOGICAL LQUERY 00070 INTEGER I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB, 00071 $ NBMIN, NN 00072 * .. 00073 * .. External Functions .. 00074 INTEGER ILAENV 00075 EXTERNAL ILAENV 00076 * .. 00077 * .. External Subroutines .. 00078 EXTERNAL DGEMM, DGEMV, DSWAP, DTRSM, DTRTRI, XERBLA 00079 * .. 00080 * .. Intrinsic Functions .. 00081 INTRINSIC MAX, MIN 00082 * .. 00083 * .. Executable Statements .. 00084 * 00085 * Test the input parameters. 00086 * 00087 INFO = 0 00088 NB = ILAENV( 1, 'DGETRI', ' ', N, -1, -1, -1 ) 00089 LWKOPT = N*NB 00090 WORK( 1 ) = LWKOPT 00091 LQUERY = ( LWORK.EQ.-1 ) 00092 IF( N.LT.0 ) THEN 00093 INFO = -1 00094 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 00095 INFO = -3 00096 ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN 00097 INFO = -6 00098 END IF 00099 IF( INFO.NE.0 ) THEN 00100 CALL XERBLA( 'DGETRI', -INFO ) 00101 RETURN 00102 ELSE IF( LQUERY ) THEN 00103 RETURN 00104 END IF 00105 * 00106 * Quick return if possible 00107 * 00108 IF( N.EQ.0 ) 00109 $ RETURN 00110 * 00111 * Form inv(U). If INFO > 0 from DTRTRI, then U is singular, 00112 * and the inverse is not computed. 00113 * 00114 CALL DTRTRI( 'Upper', 'Non-unit', N, A, LDA, INFO ) 00115 IF( INFO.GT.0 ) 00116 $ RETURN 00117 * 00118 NBMIN = 2 00119 LDWORK = N 00120 IF( NB.GT.1 .AND. NB.LT.N ) THEN 00121 IWS = MAX( LDWORK*NB, 1 ) 00122 IF( LWORK.LT.IWS ) THEN 00123 NB = LWORK / LDWORK 00124 NBMIN = MAX( 2, ILAENV( 2, 'DGETRI', ' ', N, -1, -1, -1 ) ) 00125 END IF 00126 ELSE 00127 IWS = N 00128 END IF 00129 * 00130 * Solve the equation inv(A)*L = inv(U) for inv(A). 00131 * 00132 IF( NB.LT.NBMIN .OR. NB.GE.N ) THEN 00133 * 00134 * Use unblocked code. 00135 * 00136 DO 20 J = N, 1, -1 00137 * 00138 * Copy current column of L to WORK and replace with zeros. 00139 * 00140 DO 10 I = J + 1, N 00141 WORK( I ) = A( I, J ) 00142 A( I, J ) = ZERO 00143 10 CONTINUE 00144 * 00145 * Compute current column of inv(A). 00146 * 00147 IF( J.LT.N ) 00148 $ CALL DGEMV( 'No transpose', N, N-J, -ONE, A( 1, J+1 ), 00149 $ LDA, WORK( J+1 ), 1, ONE, A( 1, J ), 1 ) 00150 20 CONTINUE 00151 ELSE 00152 * 00153 * Use blocked code. 00154 * 00155 NN = ( ( N-1 ) / NB )*NB + 1 00156 DO 50 J = NN, 1, -NB 00157 JB = MIN( NB, N-J+1 ) 00158 * 00159 * Copy current block column of L to WORK and replace with 00160 * zeros. 00161 * 00162 DO 40 JJ = J, J + JB - 1 00163 DO 30 I = JJ + 1, N 00164 WORK( I+( JJ-J )*LDWORK ) = A( I, JJ ) 00165 A( I, JJ ) = ZERO 00166 30 CONTINUE 00167 40 CONTINUE 00168 * 00169 * Compute current block column of inv(A). 00170 * 00171 IF( J+JB.LE.N ) 00172 $ CALL DGEMM( 'No transpose', 'No transpose', N, JB, 00173 $ N-J-JB+1, -ONE, A( 1, J+JB ), LDA, 00174 $ WORK( J+JB ), LDWORK, ONE, A( 1, J ), LDA ) 00175 CALL DTRSM( 'Right', 'Lower', 'No transpose', 'Unit', N, JB, 00176 $ ONE, WORK( J ), LDWORK, A( 1, J ), LDA ) 00177 50 CONTINUE 00178 END IF 00179 * 00180 * Apply column interchanges. 00181 * 00182 DO 60 J = N - 1, 1, -1 00183 JP = IPIV( J ) 00184 IF( JP.NE.J ) 00185 $ CALL DSWAP( N, A( 1, J ), 1, A( 1, JP ), 1 ) 00186 60 CONTINUE 00187 * 00188 WORK( 1 ) = IWS 00189 RETURN 00190 * 00191 * End of DGETRI 00192 * 00193 END