1*> \brief \b DGETRI 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8*> \htmlonly 9*> Download DGETRI + dependencies 10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgetri.f"> 11*> [TGZ]</a> 12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgetri.f"> 13*> [ZIP]</a> 14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgetri.f"> 15*> [TXT]</a> 16*> \endhtmlonly 17* 18* Definition: 19* =========== 20* 21* SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO ) 22* 23* .. Scalar Arguments .. 24* INTEGER INFO, LDA, LWORK, N 25* .. 26* .. Array Arguments .. 27* INTEGER IPIV( * ) 28* DOUBLE PRECISION A( LDA, * ), WORK( * ) 29* .. 30* 31* 32*> \par Purpose: 33* ============= 34*> 35*> \verbatim 36*> 37*> DGETRI computes the inverse of a matrix using the LU factorization 38*> computed by DGETRF. 39*> 40*> This method inverts U and then computes inv(A) by solving the system 41*> inv(A)*L = inv(U) for inv(A). 42*> \endverbatim 43* 44* Arguments: 45* ========== 46* 47*> \param[in] N 48*> \verbatim 49*> N is INTEGER 50*> The order of the matrix A. N >= 0. 51*> \endverbatim 52*> 53*> \param[in,out] A 54*> \verbatim 55*> A is DOUBLE PRECISION array, dimension (LDA,N) 56*> On entry, the factors L and U from the factorization 57*> A = P*L*U as computed by DGETRF. 58*> On exit, if INFO = 0, the inverse of the original matrix A. 59*> \endverbatim 60*> 61*> \param[in] LDA 62*> \verbatim 63*> LDA is INTEGER 64*> The leading dimension of the array A. LDA >= max(1,N). 65*> \endverbatim 66*> 67*> \param[in] IPIV 68*> \verbatim 69*> IPIV is INTEGER array, dimension (N) 70*> The pivot indices from DGETRF; for 1<=i<=N, row i of the 71*> matrix was interchanged with row IPIV(i). 72*> \endverbatim 73*> 74*> \param[out] WORK 75*> \verbatim 76*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) 77*> On exit, if INFO=0, then WORK(1) returns the optimal LWORK. 78*> \endverbatim 79*> 80*> \param[in] LWORK 81*> \verbatim 82*> LWORK is INTEGER 83*> The dimension of the array WORK. LWORK >= max(1,N). 84*> For optimal performance LWORK >= N*NB, where NB is 85*> the optimal blocksize returned by ILAENV. 86*> 87*> If LWORK = -1, then a workspace query is assumed; the routine 88*> only calculates the optimal size of the WORK array, returns 89*> this value as the first entry of the WORK array, and no error 90*> message related to LWORK is issued by XERBLA. 91*> \endverbatim 92*> 93*> \param[out] INFO 94*> \verbatim 95*> INFO is INTEGER 96*> = 0: successful exit 97*> < 0: if INFO = -i, the i-th argument had an illegal value 98*> > 0: if INFO = i, U(i,i) is exactly zero; the matrix is 99*> singular and its inverse could not be computed. 100*> \endverbatim 101* 102* Authors: 103* ======== 104* 105*> \author Univ. of Tennessee 106*> \author Univ. of California Berkeley 107*> \author Univ. of Colorado Denver 108*> \author NAG Ltd. 109* 110*> \ingroup doubleGEcomputational 111* 112* ===================================================================== 113 SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO ) 114* 115* -- LAPACK computational routine -- 116* -- LAPACK is a software package provided by Univ. of Tennessee, -- 117* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 118* 119* .. Scalar Arguments .. 120 INTEGER INFO, LDA, LWORK, N 121* .. 122* .. Array Arguments .. 123 INTEGER IPIV( * ) 124 DOUBLE PRECISION A( LDA, * ), WORK( * ) 125* .. 126* 127* ===================================================================== 128* 129* .. Parameters .. 130 DOUBLE PRECISION ZERO, ONE 131 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 132* .. 133* .. Local Scalars .. 134 LOGICAL LQUERY 135 INTEGER I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB, 136 $ NBMIN, NN 137* .. 138* .. External Functions .. 139 INTEGER ILAENV 140 EXTERNAL ILAENV 141* .. 142* .. External Subroutines .. 143 EXTERNAL DGEMM, DGEMV, DSWAP, DTRSM, DTRTRI, XERBLA 144* .. 145* .. Intrinsic Functions .. 146 INTRINSIC MAX, MIN 147* .. 148* .. Executable Statements .. 149* 150* Test the input parameters. 151* 152 INFO = 0 153 NB = ILAENV( 1, 'DGETRI', ' ', N, -1, -1, -1 ) 154 LWKOPT = N*NB 155 WORK( 1 ) = LWKOPT 156 LQUERY = ( LWORK.EQ.-1 ) 157 IF( N.LT.0 ) THEN 158 INFO = -1 159 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 160 INFO = -3 161 ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN 162 INFO = -6 163 END IF 164 IF( INFO.NE.0 ) THEN 165 CALL XERBLA( 'DGETRI', -INFO ) 166 RETURN 167 ELSE IF( LQUERY ) THEN 168 RETURN 169 END IF 170* 171* Quick return if possible 172* 173 IF( N.EQ.0 ) 174 $ RETURN 175* 176* Form inv(U). If INFO > 0 from DTRTRI, then U is singular, 177* and the inverse is not computed. 178* 179 CALL DTRTRI( 'Upper', 'Non-unit', N, A, LDA, INFO ) 180 IF( INFO.GT.0 ) 181 $ RETURN 182* 183 NBMIN = 2 184 LDWORK = N 185 IF( NB.GT.1 .AND. NB.LT.N ) THEN 186 IWS = MAX( LDWORK*NB, 1 ) 187 IF( LWORK.LT.IWS ) THEN 188 NB = LWORK / LDWORK 189 NBMIN = MAX( 2, ILAENV( 2, 'DGETRI', ' ', N, -1, -1, -1 ) ) 190 END IF 191 ELSE 192 IWS = N 193 END IF 194* 195* Solve the equation inv(A)*L = inv(U) for inv(A). 196* 197 IF( NB.LT.NBMIN .OR. NB.GE.N ) THEN 198* 199* Use unblocked code. 200* 201 DO 20 J = N, 1, -1 202* 203* Copy current column of L to WORK and replace with zeros. 204* 205 DO 10 I = J + 1, N 206 WORK( I ) = A( I, J ) 207 A( I, J ) = ZERO 208 10 CONTINUE 209* 210* Compute current column of inv(A). 211* 212 IF( J.LT.N ) 213 $ CALL DGEMV( 'No transpose', N, N-J, -ONE, A( 1, J+1 ), 214 $ LDA, WORK( J+1 ), 1, ONE, A( 1, J ), 1 ) 215 20 CONTINUE 216 ELSE 217* 218* Use blocked code. 219* 220 NN = ( ( N-1 ) / NB )*NB + 1 221 DO 50 J = NN, 1, -NB 222 JB = MIN( NB, N-J+1 ) 223* 224* Copy current block column of L to WORK and replace with 225* zeros. 226* 227 DO 40 JJ = J, J + JB - 1 228 DO 30 I = JJ + 1, N 229 WORK( I+( JJ-J )*LDWORK ) = A( I, JJ ) 230 A( I, JJ ) = ZERO 231 30 CONTINUE 232 40 CONTINUE 233* 234* Compute current block column of inv(A). 235* 236 IF( J+JB.LE.N ) 237 $ CALL DGEMM( 'No transpose', 'No transpose', N, JB, 238 $ N-J-JB+1, -ONE, A( 1, J+JB ), LDA, 239 $ WORK( J+JB ), LDWORK, ONE, A( 1, J ), LDA ) 240 CALL DTRSM( 'Right', 'Lower', 'No transpose', 'Unit', N, JB, 241 $ ONE, WORK( J ), LDWORK, A( 1, J ), LDA ) 242 50 CONTINUE 243 END IF 244* 245* Apply column interchanges. 246* 247 DO 60 J = N - 1, 1, -1 248 JP = IPIV( J ) 249 IF( JP.NE.J ) 250 $ CALL DSWAP( N, A( 1, J ), 1, A( 1, JP ), 1 ) 251 60 CONTINUE 252* 253 WORK( 1 ) = IWS 254 RETURN 255* 256* End of DGETRI 257* 258 END 259