LAPACK 3.12.0
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 113 of file zgetri.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 COMPLEX*16 A( LDA, * ), WORK( * )
125* ..
126*
127* =====================================================================
128*
129* .. Parameters ..
130 COMPLEX*16 ZERO, ONE
131 parameter( zero = ( 0.0d+0, 0.0d+0 ),
132 $ one = ( 1.0d+0, 0.0d+0 ) )
133* ..
134* .. Local Scalars ..
135 LOGICAL LQUERY
136 INTEGER I, IWS, J, JB, JJ, JP, LDWORK, LWKOPT, NB,
137 $ NBMIN, NN
138* ..
139* .. External Functions ..
140 INTEGER ILAENV
141 EXTERNAL ilaenv
142* ..
143* .. External Subroutines ..
144 EXTERNAL xerbla, zgemm, zgemv, zswap, ztrsm, ztrtri
145* ..
146* .. Intrinsic Functions ..
147 INTRINSIC max, min
148* ..
149* .. Executable Statements ..
150*
151* Test the input parameters.
152*
153 info = 0
154 nb = ilaenv( 1, 'ZGETRI', ' ', n, -1, -1, -1 )
155 lwkopt = n*nb
156 work( 1 ) = lwkopt
157 lquery = ( lwork.EQ.-1 )
158 IF( n.LT.0 ) THEN
159 info = -1
160 ELSE IF( lda.LT.max( 1, n ) ) THEN
161 info = -3
162 ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
163 info = -6
164 END IF
165 IF( info.NE.0 ) THEN
166 CALL xerbla( 'ZGETRI', -info )
167 RETURN
168 ELSE IF( lquery ) THEN
169 RETURN
170 END IF
171*
172* Quick return if possible
173*
174 IF( n.EQ.0 )
175 $ RETURN
176*
177* Form inv(U). If INFO > 0 from ZTRTRI, then U is singular,
178* and the inverse is not computed.
179*
180 CALL ztrtri( 'Upper', 'Non-unit', n, a, lda, info )
181 IF( info.GT.0 )
182 $ RETURN
183*
184 nbmin = 2
185 ldwork = n
186 IF( nb.GT.1 .AND. nb.LT.n ) THEN
187 iws = max( ldwork*nb, 1 )
188 IF( lwork.LT.iws ) THEN
189 nb = lwork / ldwork
190 nbmin = max( 2, ilaenv( 2, 'ZGETRI', ' ', n, -1, -1, -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, jb,
242 $ one, work( j ), ldwork, a( 1, j ), lda )
243 50 CONTINUE
244 END IF
245*
246* Apply column interchanges.
247*
248 DO 60 j = n - 1, 1, -1
249 jp = ipiv( j )
250 IF( jp.NE.j )
251 $ CALL zswap( n, a( 1, j ), 1, a( 1, jp ), 1 )
252 60 CONTINUE
253*
254 work( 1 ) = iws
255 RETURN
256*
257* End of ZGETRI
258*
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:162
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:109
Here is the call graph for this function:
Here is the caller graph for this function: