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

◆ zgetri()

subroutine zgetri ( integer n,
complex*16, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) ipiv,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZGETRI

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

Purpose:
!>
!> ZGETRI computes the inverse of a matrix using the LU factorization
!> computed by ZGETRF.
!>
!> 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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the factors L and U from the factorization
!>          A = P*L*U as computed by ZGETRF.
!>          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 ZGETRF; for 1<=i<=N, row i of the
!>          matrix was interchanged with row IPIV(i).
!> 
[out]WORK
!>          WORK is COMPLEX*16 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 111 of file zgetri.f.

112*
113* -- LAPACK computational routine --
114* -- LAPACK is a software package provided by Univ. of Tennessee, --
115* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
116*
117* .. Scalar Arguments ..
118 INTEGER INFO, LDA, LWORK, N
119* ..
120* .. Array Arguments ..
121 INTEGER IPIV( * )
122 COMPLEX*16 A( LDA, * ), WORK( * )
123* ..
124*
125* =====================================================================
126*
127* .. Parameters ..
128 COMPLEX*16 ZERO, ONE
129 parameter( zero = ( 0.0d+0, 0.0d+0 ),
130 $ one = ( 1.0d+0, 0.0d+0 ) )
131* ..
132* .. Local Scalars ..
133 LOGICAL LQUERY
134 INTEGER I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB,
135 $ NBMIN, NN
136* ..
137* .. External Functions ..
138 INTEGER ILAENV
139 EXTERNAL ilaenv
140* ..
141* .. External Subroutines ..
142 EXTERNAL xerbla, zgemm, zgemv, zswap, ztrsm,
143 $ ztrtri
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, 'ZGETRI', ' ', n, -1, -1, -1 )
154 lwkopt = max( 1, 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( 'ZGETRI', -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 ZTRTRI, then U is singular,
177* and the inverse is not computed.
178*
179 CALL ztrtri( '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, 'ZGETRI', ' ', n, -1, -1,
190 $ -1 ) )
191 END IF
192 ELSE
193 iws = n
194 END IF
195*
196* Solve the equation inv(A)*L = inv(U) for inv(A).
197*
198 IF( nb.LT.nbmin .OR. nb.GE.n ) THEN
199*
200* Use unblocked code.
201*
202 DO 20 j = n, 1, -1
203*
204* Copy current column of L to WORK and replace with zeros.
205*
206 DO 10 i = j + 1, n
207 work( i ) = a( i, j )
208 a( i, j ) = zero
209 10 CONTINUE
210*
211* Compute current column of inv(A).
212*
213 IF( j.LT.n )
214 $ CALL zgemv( 'No transpose', n, n-j, -one, a( 1, j+1 ),
215 $ lda, work( j+1 ), 1, one, a( 1, j ), 1 )
216 20 CONTINUE
217 ELSE
218*
219* Use blocked code.
220*
221 nn = ( ( n-1 ) / nb )*nb + 1
222 DO 50 j = nn, 1, -nb
223 jb = min( nb, n-j+1 )
224*
225* Copy current block column of L to WORK and replace with
226* zeros.
227*
228 DO 40 jj = j, j + jb - 1
229 DO 30 i = jj + 1, n
230 work( i+( jj-j )*ldwork ) = a( i, jj )
231 a( i, jj ) = zero
232 30 CONTINUE
233 40 CONTINUE
234*
235* Compute current block column of inv(A).
236*
237 IF( j+jb.LE.n )
238 $ CALL zgemm( 'No transpose', 'No transpose', n, jb,
239 $ n-j-jb+1, -one, a( 1, j+jb ), lda,
240 $ work( j+jb ), ldwork, one, a( 1, j ), lda )
241 CALL ztrsm( 'Right', 'Lower', 'No transpose', 'Unit', n,
242 $ jb,
243 $ one, work( j ), ldwork, a( 1, j ), lda )
244 50 CONTINUE
245 END IF
246*
247* Apply column interchanges.
248*
249 DO 60 j = n - 1, 1, -1
250 jp = ipiv( j )
251 IF( jp.NE.j )
252 $ CALL zswap( n, a( 1, j ), 1, a( 1, jp ), 1 )
253 60 CONTINUE
254*
255 work( 1 ) = iws
256 RETURN
257*
258* End of ZGETRI
259*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180
subroutine ztrtri(uplo, diag, n, a, lda, info)
ZTRTRI
Definition ztrtri.f:107
Here is the call graph for this function:
Here is the caller graph for this function: