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

◆ cgetri()

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

CGETRI

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

Purpose:
 CGETRI computes the inverse of a matrix using the LU factorization
 computed by CGETRF.

 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 array, dimension (LDA,N)
          On entry, the factors L and U from the factorization
          A = P*L*U as computed by CGETRF.
          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 CGETRF; for 1<=i<=N, row i of the
          matrix was interchanged with row IPIV(i).
[out]WORK
          WORK is COMPLEX 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 cgetri.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 A( LDA, * ), WORK( * )
125* ..
126*
127* =====================================================================
128*
129* .. Parameters ..
130 COMPLEX ZERO, ONE
131 parameter( zero = ( 0.0e+0, 0.0e+0 ),
132 $ one = ( 1.0e+0, 0.0e+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 REAL SROUNDUP_LWORK
142 EXTERNAL ilaenv, sroundup_lwork
143* ..
144* .. External Subroutines ..
145 EXTERNAL cgemm, cgemv, cswap, ctrsm, ctrtri, xerbla
146* ..
147* .. Intrinsic Functions ..
148 INTRINSIC max, min
149* ..
150* .. Executable Statements ..
151*
152* Test the input parameters.
153*
154 info = 0
155 nb = ilaenv( 1, 'CGETRI', ' ', n, -1, -1, -1 )
156 lwkopt = n*nb
157 work( 1 ) = sroundup_lwork(lwkopt)
158 lquery = ( lwork.EQ.-1 )
159 IF( n.LT.0 ) THEN
160 info = -1
161 ELSE IF( lda.LT.max( 1, n ) ) THEN
162 info = -3
163 ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
164 info = -6
165 END IF
166 IF( info.NE.0 ) THEN
167 CALL xerbla( 'CGETRI', -info )
168 RETURN
169 ELSE IF( lquery ) THEN
170 RETURN
171 END IF
172*
173* Quick return if possible
174*
175 IF( n.EQ.0 )
176 $ RETURN
177*
178* Form inv(U). If INFO > 0 from CTRTRI, then U is singular,
179* and the inverse is not computed.
180*
181 CALL ctrtri( 'Upper', 'Non-unit', n, a, lda, info )
182 IF( info.GT.0 )
183 $ RETURN
184*
185 nbmin = 2
186 ldwork = n
187 IF( nb.GT.1 .AND. nb.LT.n ) THEN
188 iws = max( ldwork*nb, 1 )
189 IF( lwork.LT.iws ) THEN
190 nb = lwork / ldwork
191 nbmin = max( 2, ilaenv( 2, 'CGETRI', ' ', n, -1, -1, -1 ) )
192 END IF
193 ELSE
194 iws = n
195 END IF
196*
197* Solve the equation inv(A)*L = inv(U) for inv(A).
198*
199 IF( nb.LT.nbmin .OR. nb.GE.n ) THEN
200*
201* Use unblocked code.
202*
203 DO 20 j = n, 1, -1
204*
205* Copy current column of L to WORK and replace with zeros.
206*
207 DO 10 i = j + 1, n
208 work( i ) = a( i, j )
209 a( i, j ) = zero
210 10 CONTINUE
211*
212* Compute current column of inv(A).
213*
214 IF( j.LT.n )
215 $ CALL cgemv( 'No transpose', n, n-j, -one, a( 1, j+1 ),
216 $ lda, work( j+1 ), 1, one, a( 1, j ), 1 )
217 20 CONTINUE
218 ELSE
219*
220* Use blocked code.
221*
222 nn = ( ( n-1 ) / nb )*nb + 1
223 DO 50 j = nn, 1, -nb
224 jb = min( nb, n-j+1 )
225*
226* Copy current block column of L to WORK and replace with
227* zeros.
228*
229 DO 40 jj = j, j + jb - 1
230 DO 30 i = jj + 1, n
231 work( i+( jj-j )*ldwork ) = a( i, jj )
232 a( i, jj ) = zero
233 30 CONTINUE
234 40 CONTINUE
235*
236* Compute current block column of inv(A).
237*
238 IF( j+jb.LE.n )
239 $ CALL cgemm( 'No transpose', 'No transpose', n, jb,
240 $ n-j-jb+1, -one, a( 1, j+jb ), lda,
241 $ work( j+jb ), ldwork, one, a( 1, j ), lda )
242 CALL ctrsm( 'Right', 'Lower', 'No transpose', 'Unit', n, 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 cswap( n, a( 1, j ), 1, a( 1, jp ), 1 )
253 60 CONTINUE
254*
255 work( 1 ) = sroundup_lwork(iws)
256 RETURN
257*
258* End of CGETRI
259*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:160
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
subroutine ctrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRSM
Definition ctrsm.f:180
subroutine ctrtri(uplo, diag, n, a, lda, info)
CTRTRI
Definition ctrtri.f:109
Here is the call graph for this function:
Here is the caller graph for this function: