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

◆ zgbcon()

subroutine zgbcon ( character norm,
integer n,
integer kl,
integer ku,
complex*16, dimension( ldab, * ) ab,
integer ldab,
integer, dimension( * ) ipiv,
double precision anorm,
double precision rcond,
complex*16, dimension( * ) work,
double precision, dimension( * ) rwork,
integer info )

ZGBCON

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

Purpose:
!>
!> ZGBCON estimates the reciprocal of the condition number of a complex
!> general band matrix A, in either the 1-norm or the infinity-norm,
!> using the LU factorization computed by ZGBTRF.
!>
!> An estimate is obtained for norm(inv(A)), and the reciprocal of the
!> condition number is computed as
!>    RCOND = 1 / ( norm(A) * norm(inv(A)) ).
!> 
Parameters
[in]NORM
!>          NORM is CHARACTER*1
!>          Specifies whether the 1-norm condition number or the
!>          infinity-norm condition number is required:
!>          = '1' or 'O':  1-norm;
!>          = 'I':         Infinity-norm.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in]KL
!>          KL is INTEGER
!>          The number of subdiagonals within the band of A.  KL >= 0.
!> 
[in]KU
!>          KU is INTEGER
!>          The number of superdiagonals within the band of A.  KU >= 0.
!> 
[in]AB
!>          AB is COMPLEX*16 array, dimension (LDAB,N)
!>          Details of the LU factorization of the band matrix A, as
!>          computed by ZGBTRF.  U is stored as an upper triangular band
!>          matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and
!>          the multipliers used during the factorization are stored in
!>          rows KL+KU+2 to 2*KL+KU+1.
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= 2*KL+KU+1.
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          The pivot indices; for 1 <= i <= N, row i of the matrix was
!>          interchanged with row IPIV(i).
!> 
[in]ANORM
!>          ANORM is DOUBLE PRECISION
!>          If NORM = '1' or 'O', the 1-norm of the original matrix A.
!>          If NORM = 'I', the infinity-norm of the original matrix A.
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>          The reciprocal of the condition number of the matrix A,
!>          computed as RCOND = 1/(norm(A) * norm(inv(A))).
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (2*N)
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 143 of file zgbcon.f.

146*
147* -- LAPACK computational routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 CHARACTER NORM
153 INTEGER INFO, KL, KU, LDAB, N
154 DOUBLE PRECISION ANORM, RCOND
155* ..
156* .. Array Arguments ..
157 INTEGER IPIV( * )
158 DOUBLE PRECISION RWORK( * )
159 COMPLEX*16 AB( LDAB, * ), WORK( * )
160* ..
161*
162* =====================================================================
163*
164* .. Parameters ..
165 DOUBLE PRECISION ONE, ZERO
166 parameter( one = 1.0d+0, zero = 0.0d+0 )
167* ..
168* .. Local Scalars ..
169 LOGICAL LNOTI, ONENRM
170 CHARACTER NORMIN
171 INTEGER IX, J, JP, KASE, KASE1, KD, LM
172 DOUBLE PRECISION AINVNM, SCALE, SMLNUM
173 COMPLEX*16 T, ZDUM
174* ..
175* .. Local Arrays ..
176 INTEGER ISAVE( 3 )
177* ..
178* .. External Functions ..
179 LOGICAL LSAME
180 INTEGER IZAMAX
181 DOUBLE PRECISION DLAMCH
182 COMPLEX*16 ZDOTC
183 EXTERNAL lsame, izamax, dlamch, zdotc
184* ..
185* .. External Subroutines ..
186 EXTERNAL xerbla, zaxpy, zdrscl, zlacn2,
187 $ zlatbs
188* ..
189* .. Intrinsic Functions ..
190 INTRINSIC abs, dble, dimag, min
191* ..
192* .. Statement Functions ..
193 DOUBLE PRECISION CABS1
194* ..
195* .. Statement Function definitions ..
196 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
197* ..
198* .. Executable Statements ..
199*
200* Test the input parameters.
201*
202 info = 0
203 onenrm = norm.EQ.'1' .OR. lsame( norm, 'O' )
204 IF( .NOT.onenrm .AND. .NOT.lsame( norm, 'I' ) ) THEN
205 info = -1
206 ELSE IF( n.LT.0 ) THEN
207 info = -2
208 ELSE IF( kl.LT.0 ) THEN
209 info = -3
210 ELSE IF( ku.LT.0 ) THEN
211 info = -4
212 ELSE IF( ldab.LT.2*kl+ku+1 ) THEN
213 info = -6
214 ELSE IF( anorm.LT.zero ) THEN
215 info = -8
216 END IF
217 IF( info.NE.0 ) THEN
218 CALL xerbla( 'ZGBCON', -info )
219 RETURN
220 END IF
221*
222* Quick return if possible
223*
224 rcond = zero
225 IF( n.EQ.0 ) THEN
226 rcond = one
227 RETURN
228 ELSE IF( anorm.EQ.zero ) THEN
229 RETURN
230 END IF
231*
232 smlnum = dlamch( 'Safe minimum' )
233*
234* Estimate the norm of inv(A).
235*
236 ainvnm = zero
237 normin = 'N'
238 IF( onenrm ) THEN
239 kase1 = 1
240 ELSE
241 kase1 = 2
242 END IF
243 kd = kl + ku + 1
244 lnoti = kl.GT.0
245 kase = 0
246 10 CONTINUE
247 CALL zlacn2( n, work( n+1 ), work, ainvnm, kase, isave )
248 IF( kase.NE.0 ) THEN
249 IF( kase.EQ.kase1 ) THEN
250*
251* Multiply by inv(L).
252*
253 IF( lnoti ) THEN
254 DO 20 j = 1, n - 1
255 lm = min( kl, n-j )
256 jp = ipiv( j )
257 t = work( jp )
258 IF( jp.NE.j ) THEN
259 work( jp ) = work( j )
260 work( j ) = t
261 END IF
262 CALL zaxpy( lm, -t, ab( kd+1, j ), 1, work( j+1 ),
263 $ 1 )
264 20 CONTINUE
265 END IF
266*
267* Multiply by inv(U).
268*
269 CALL zlatbs( 'Upper', 'No transpose', 'Non-unit', normin,
270 $ n,
271 $ kl+ku, ab, ldab, work, scale, rwork, info )
272 ELSE
273*
274* Multiply by inv(U**H).
275*
276 CALL zlatbs( 'Upper', 'Conjugate transpose', 'Non-unit',
277 $ normin, n, kl+ku, ab, ldab, work, scale, rwork,
278 $ info )
279*
280* Multiply by inv(L**H).
281*
282 IF( lnoti ) THEN
283 DO 30 j = n - 1, 1, -1
284 lm = min( kl, n-j )
285 work( j ) = work( j ) - zdotc( lm, ab( kd+1, j ),
286 $ 1,
287 $ work( j+1 ), 1 )
288 jp = ipiv( j )
289 IF( jp.NE.j ) THEN
290 t = work( jp )
291 work( jp ) = work( j )
292 work( j ) = t
293 END IF
294 30 CONTINUE
295 END IF
296 END IF
297*
298* Divide X by 1/SCALE if doing so will not cause overflow.
299*
300 normin = 'Y'
301 IF( scale.NE.one ) THEN
302 ix = izamax( n, work, 1 )
303 IF( scale.LT.cabs1( work( ix ) )*smlnum .OR. scale.EQ.zero )
304 $ GO TO 40
305 CALL zdrscl( n, scale, work, 1 )
306 END IF
307 GO TO 10
308 END IF
309*
310* Compute the estimate of the reciprocal condition number.
311*
312 IF( ainvnm.NE.zero )
313 $ rcond = ( one / ainvnm ) / anorm
314*
315 40 CONTINUE
316 RETURN
317*
318* End of ZGBCON
319*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
complex *16 function zdotc(n, zx, incx, zy, incy)
ZDOTC
Definition zdotc.f:83
integer function izamax(n, zx, incx)
IZAMAX
Definition izamax.f:71
subroutine zlacn2(n, v, x, est, kase, isave)
ZLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition zlacn2.f:131
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine zlatbs(uplo, trans, diag, normin, n, kd, ab, ldab, x, scale, cnorm, info)
ZLATBS solves a triangular banded system of equations.
Definition zlatbs.f:242
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zdrscl(n, sa, sx, incx)
ZDRSCL multiplies a vector by the reciprocal of a real scalar.
Definition zdrscl.f:82
Here is the call graph for this function:
Here is the caller graph for this function: