LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dgetri()

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.

Definition at line 113 of file dgetri.f.

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*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
subroutine dtrtri(uplo, diag, n, a, lda, info)
DTRTRI
Definition dtrtri.f:109
Here is the call graph for this function:
Here is the caller graph for this function: