3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DGETRI + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgetri.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgetri.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgetri.f">
21 * SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO )
23 * .. Scalar Arguments ..
24 * INTEGER INFO, LDA, LWORK, N
26 * .. Array Arguments ..
28 * DOUBLE PRECISION A( LDA, * ), WORK( * )
37 *> DGETRI computes the inverse of a matrix using the LU factorization
38 *> computed by DGETRF.
40 *> This method inverts U and then computes inv(A) by solving the system
41 *> inv(A)*L = inv(U) for inv(A).
50 *> The order of the matrix A. N >= 0.
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.
64 *> The leading dimension of the array A. LDA >= max(1,N).
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).
76 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
77 *> On exit, if INFO=0, then WORK(1) returns the optimal LWORK.
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.
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.
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.
105 *> \author Univ. of Tennessee
106 *> \author Univ. of California Berkeley
107 *> \author Univ. of Colorado Denver
110 *> \date November 2011
112 *> \ingroup doubleGEcomputational
114 * =====================================================================
115 SUBROUTINE DGETRI( N, A, LDA, IPIV, WORK, LWORK, INFO )
117 * -- LAPACK computational routine (version 3.4.0) --
118 * -- LAPACK is a software package provided by Univ. of Tennessee, --
119 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
122 * .. Scalar Arguments ..
123 INTEGER INFO, LDA, LWORK, N
125 * .. Array Arguments ..
127 DOUBLE PRECISION A( LDA, * ), WORK( * )
130 * =====================================================================
133 DOUBLE PRECISION ZERO, ONE
134 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
136 * .. Local Scalars ..
138 INTEGER I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB,
141 * .. External Functions ..
145 * .. External Subroutines ..
146 EXTERNAL DGEMM, DGEMV, DSWAP, DTRSM, DTRTRI, XERBLA
148 * .. Intrinsic Functions ..
151 * .. Executable Statements ..
153 * Test the input parameters.
156 NB = ILAENV( 1, 'DGETRI', ' ', N, -1, -1, -1 )
159 LQUERY = ( LWORK.EQ.-1 )
162 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
164 ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN
168 CALL XERBLA( 'DGETRI', -INFO )
170 ELSE IF( LQUERY ) THEN
174 * Quick return if possible
179 * Form inv(U). If INFO > 0 from DTRTRI, then U is singular,
180 * and the inverse is not computed.
182 CALL DTRTRI( 'Upper', 'Non-unit', N, A, LDA, INFO )
188 IF( NB.GT.1 .AND. NB.LT.N ) THEN
189 IWS = MAX( LDWORK*NB, 1 )
190 IF( LWORK.LT.IWS ) THEN
192 NBMIN = MAX( 2, ILAENV( 2, 'DGETRI', ' ', N, -1, -1, -1 ) )
198 * Solve the equation inv(A)*L = inv(U) for inv(A).
200 IF( NB.LT.NBMIN .OR. NB.GE.N ) THEN
202 * Use unblocked code.
206 * Copy current column of L to WORK and replace with zeros.
209 WORK( I ) = A( I, J )
213 * Compute current column of inv(A).
216 $ CALL DGEMV( 'No transpose', N, N-J, -ONE, A( 1, J+1 ),
217 $ LDA, WORK( J+1 ), 1, ONE, A( 1, J ), 1 )
223 NN = ( ( N-1 ) / NB )*NB + 1
225 JB = MIN( NB, N-J+1 )
227 * Copy current block column of L to WORK and replace with
230 DO 40 JJ = J, J + JB - 1
232 WORK( I+( JJ-J )*LDWORK ) = A( I, JJ )
237 * Compute current block column of inv(A).
240 $ CALL DGEMM( 'No transpose', 'No transpose', N, JB,
241 $ N-J-JB+1, -ONE, A( 1, J+JB ), LDA,
242 $ WORK( J+JB ), LDWORK, ONE, A( 1, J ), LDA )
243 CALL DTRSM( 'Right', 'Lower', 'No transpose', 'Unit', N, JB,
244 $ ONE, WORK( J ), LDWORK, A( 1, J ), LDA )
248 * Apply column interchanges.
250 DO 60 J = N - 1, 1, -1
253 $ CALL DSWAP( N, A( 1, J ), 1, A( 1, JP ), 1 )