LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgetri ( integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
double precision, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

DGETRI

Download DGETRI + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 DGETRI computes the inverse of a matrix using the LU factorization
 computed by DGETRF.

 This method inverts U and then computes inv(A) by solving the system
 inv(A)*L = inv(U) for inv(A).
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the factors L and U from the factorization
          A = P*L*U as computed by DGETRF.
          On exit, if INFO = 0, the inverse of the original matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          The pivot indices from DGETRF; for 1<=i<=N, row i of the
          matrix was interchanged with row IPIV(i).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
          On exit, if INFO=0, then WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= max(1,N).
          For optimal performance LWORK >= N*NB, where NB is
          the optimal blocksize returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, U(i,i) is exactly zero; the matrix is
                singular and its inverse could not be computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 116 of file dgetri.f.

116 *
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..--
120 * November 2011
121 *
122 * .. Scalar Arguments ..
123  INTEGER info, lda, lwork, n
124 * ..
125 * .. Array Arguments ..
126  INTEGER ipiv( * )
127  DOUBLE PRECISION a( lda, * ), work( * )
128 * ..
129 *
130 * =====================================================================
131 *
132 * .. Parameters ..
133  DOUBLE PRECISION zero, one
134  parameter ( zero = 0.0d+0, one = 1.0d+0 )
135 * ..
136 * .. Local Scalars ..
137  LOGICAL lquery
138  INTEGER i, iws, j, jb, jj, jp, ldwork, lwkopt, nb,
139  $ nbmin, nn
140 * ..
141 * .. External Functions ..
142  INTEGER ilaenv
143  EXTERNAL ilaenv
144 * ..
145 * .. External Subroutines ..
146  EXTERNAL dgemm, dgemv, dswap, dtrsm, dtrtri, xerbla
147 * ..
148 * .. Intrinsic Functions ..
149  INTRINSIC max, min
150 * ..
151 * .. Executable Statements ..
152 *
153 * Test the input parameters.
154 *
155  info = 0
156  nb = ilaenv( 1, 'DGETRI', ' ', n, -1, -1, -1 )
157  lwkopt = n*nb
158  work( 1 ) = lwkopt
159  lquery = ( lwork.EQ.-1 )
160  IF( n.LT.0 ) THEN
161  info = -1
162  ELSE IF( lda.LT.max( 1, n ) ) THEN
163  info = -3
164  ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
165  info = -6
166  END IF
167  IF( info.NE.0 ) THEN
168  CALL xerbla( 'DGETRI', -info )
169  RETURN
170  ELSE IF( lquery ) THEN
171  RETURN
172  END IF
173 *
174 * Quick return if possible
175 *
176  IF( n.EQ.0 )
177  $ RETURN
178 *
179 * Form inv(U). If INFO > 0 from DTRTRI, then U is singular,
180 * and the inverse is not computed.
181 *
182  CALL dtrtri( 'Upper', 'Non-unit', n, a, lda, info )
183  IF( info.GT.0 )
184  $ RETURN
185 *
186  nbmin = 2
187  ldwork = n
188  IF( nb.GT.1 .AND. nb.LT.n ) THEN
189  iws = max( ldwork*nb, 1 )
190  IF( lwork.LT.iws ) THEN
191  nb = lwork / ldwork
192  nbmin = max( 2, ilaenv( 2, 'DGETRI', ' ', n, -1, -1, -1 ) )
193  END IF
194  ELSE
195  iws = n
196  END IF
197 *
198 * Solve the equation inv(A)*L = inv(U) for inv(A).
199 *
200  IF( nb.LT.nbmin .OR. nb.GE.n ) THEN
201 *
202 * Use unblocked code.
203 *
204  DO 20 j = n, 1, -1
205 *
206 * Copy current column of L to WORK and replace with zeros.
207 *
208  DO 10 i = j + 1, n
209  work( i ) = a( i, j )
210  a( i, j ) = zero
211  10 CONTINUE
212 *
213 * Compute current column of inv(A).
214 *
215  IF( j.LT.n )
216  $ CALL dgemv( 'No transpose', n, n-j, -one, a( 1, j+1 ),
217  $ lda, work( j+1 ), 1, one, a( 1, j ), 1 )
218  20 CONTINUE
219  ELSE
220 *
221 * Use blocked code.
222 *
223  nn = ( ( n-1 ) / nb )*nb + 1
224  DO 50 j = nn, 1, -nb
225  jb = min( nb, n-j+1 )
226 *
227 * Copy current block column of L to WORK and replace with
228 * zeros.
229 *
230  DO 40 jj = j, j + jb - 1
231  DO 30 i = jj + 1, n
232  work( i+( jj-j )*ldwork ) = a( i, jj )
233  a( i, jj ) = zero
234  30 CONTINUE
235  40 CONTINUE
236 *
237 * Compute current block column of inv(A).
238 *
239  IF( j+jb.LE.n )
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 )
245  50 CONTINUE
246  END IF
247 *
248 * Apply column interchanges.
249 *
250  DO 60 j = n - 1, 1, -1
251  jp = ipiv( j )
252  IF( jp.NE.j )
253  $ CALL dswap( n, a( 1, j ), 1, a( 1, jp ), 1 )
254  60 CONTINUE
255 *
256  work( 1 ) = iws
257  RETURN
258 *
259 * End of DGETRI
260 *
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
Definition: dtrsm.f:183
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine dtrtri(UPLO, DIAG, N, A, LDA, INFO)
DTRTRI
Definition: dtrtri.f:111

Here is the call graph for this function:

Here is the caller graph for this function: