LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
double precision function zla_gbrcond_c ( character  TRANS,
integer  N,
integer  KL,
integer  KU,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
complex*16, dimension( ldafb, * )  AFB,
integer  LDAFB,
integer, dimension( * )  IPIV,
double precision, dimension( * )  C,
logical  CAPPLY,
integer  INFO,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK 
)

ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general banded matrices.

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

Purpose:
    ZLA_GBRCOND_C Computes the infinity norm condition number of
    op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector.
Parameters
[in]TRANS
          TRANS is CHARACTER*1
     Specifies the form of the system of equations:
       = 'N':  A * X = B     (No transpose)
       = 'T':  A**T * X = B  (Transpose)
       = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)
[in]N
          N is INTEGER
     The number of linear equations, i.e., 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)
     On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
     The j-th column of A is stored in the j-th column of the
     array AB as follows:
     AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
[in]LDAB
          LDAB is INTEGER
     The leading dimension of the array AB.  LDAB >= KL+KU+1.
[in]AFB
          AFB is COMPLEX*16 array, dimension (LDAFB,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]LDAFB
          LDAFB is INTEGER
     The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
[in]IPIV
          IPIV is INTEGER array, dimension (N)
     The pivot indices from the factorization A = P*L*U
     as computed by ZGBTRF; row i of the matrix was interchanged
     with row IPIV(i).
[in]C
          C is DOUBLE PRECISION array, dimension (N)
     The vector C in the formula op(A) * inv(diag(C)).
[in]CAPPLY
          CAPPLY is LOGICAL
     If .TRUE. then access the vector C in the formula above.
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[in]WORK
          WORK is COMPLEX*16 array, dimension (2*N).
     Workspace.
[in]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 165 of file zla_gbrcond_c.f.

165 *
166 * -- LAPACK computational routine (version 3.4.2) --
167 * -- LAPACK is a software package provided by Univ. of Tennessee, --
168 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
169 * September 2012
170 *
171 * .. Scalar Arguments ..
172  CHARACTER trans
173  LOGICAL capply
174  INTEGER n, kl, ku, kd, ke, ldab, ldafb, info
175 * ..
176 * .. Array Arguments ..
177  INTEGER ipiv( * )
178  COMPLEX*16 ab( ldab, * ), afb( ldafb, * ), work( * )
179  DOUBLE PRECISION c( * ), rwork( * )
180 *
181 *
182 * =====================================================================
183 *
184 * .. Local Scalars ..
185  LOGICAL notrans
186  INTEGER kase, i, j
187  DOUBLE PRECISION ainvnm, anorm, tmp
188  COMPLEX*16 zdum
189 * ..
190 * .. Local Arrays ..
191  INTEGER isave( 3 )
192 * ..
193 * .. External Functions ..
194  LOGICAL lsame
195  EXTERNAL lsame
196 * ..
197 * .. External Subroutines ..
198  EXTERNAL zlacn2, zgbtrs, xerbla
199 * ..
200 * .. Intrinsic Functions ..
201  INTRINSIC abs, max
202 * ..
203 * .. Statement Functions ..
204  DOUBLE PRECISION cabs1
205 * ..
206 * .. Statement Function Definitions ..
207  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
208 * ..
209 * .. Executable Statements ..
210  zla_gbrcond_c = 0.0d+0
211 *
212  info = 0
213  notrans = lsame( trans, 'N' )
214  IF ( .NOT. notrans .AND. .NOT. lsame( trans, 'T' ) .AND. .NOT.
215  $ lsame( trans, 'C' ) ) THEN
216  info = -1
217  ELSE IF( n.LT.0 ) THEN
218  info = -2
219  ELSE IF( kl.LT.0 .OR. kl.GT.n-1 ) THEN
220  info = -3
221  ELSE IF( ku.LT.0 .OR. ku.GT.n-1 ) THEN
222  info = -4
223  ELSE IF( ldab.LT.kl+ku+1 ) THEN
224  info = -6
225  ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
226  info = -8
227  END IF
228  IF( info.NE.0 ) THEN
229  CALL xerbla( 'ZLA_GBRCOND_C', -info )
230  RETURN
231  END IF
232 *
233 * Compute norm of op(A)*op2(C).
234 *
235  anorm = 0.0d+0
236  kd = ku + 1
237  ke = kl + 1
238  IF ( notrans ) THEN
239  DO i = 1, n
240  tmp = 0.0d+0
241  IF ( capply ) THEN
242  DO j = max( i-kl, 1 ), min( i+ku, n )
243  tmp = tmp + cabs1( ab( kd+i-j, j ) ) / c( j )
244  END DO
245  ELSE
246  DO j = max( i-kl, 1 ), min( i+ku, n )
247  tmp = tmp + cabs1( ab( kd+i-j, j ) )
248  END DO
249  END IF
250  rwork( i ) = tmp
251  anorm = max( anorm, tmp )
252  END DO
253  ELSE
254  DO i = 1, n
255  tmp = 0.0d+0
256  IF ( capply ) THEN
257  DO j = max( i-kl, 1 ), min( i+ku, n )
258  tmp = tmp + cabs1( ab( ke-i+j, i ) ) / c( j )
259  END DO
260  ELSE
261  DO j = max( i-kl, 1 ), min( i+ku, n )
262  tmp = tmp + cabs1( ab( ke-i+j, i ) )
263  END DO
264  END IF
265  rwork( i ) = tmp
266  anorm = max( anorm, tmp )
267  END DO
268  END IF
269 *
270 * Quick return if possible.
271 *
272  IF( n.EQ.0 ) THEN
273  zla_gbrcond_c = 1.0d+0
274  RETURN
275  ELSE IF( anorm .EQ. 0.0d+0 ) THEN
276  RETURN
277  END IF
278 *
279 * Estimate the norm of inv(op(A)).
280 *
281  ainvnm = 0.0d+0
282 *
283  kase = 0
284  10 CONTINUE
285  CALL zlacn2( n, work( n+1 ), work, ainvnm, kase, isave )
286  IF( kase.NE.0 ) THEN
287  IF( kase.EQ.2 ) THEN
288 *
289 * Multiply by R.
290 *
291  DO i = 1, n
292  work( i ) = work( i ) * rwork( i )
293  END DO
294 *
295  IF ( notrans ) THEN
296  CALL zgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
297  $ ipiv, work, n, info )
298  ELSE
299  CALL zgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
300  $ ldafb, ipiv, work, n, info )
301  ENDIF
302 *
303 * Multiply by inv(C).
304 *
305  IF ( capply ) THEN
306  DO i = 1, n
307  work( i ) = work( i ) * c( i )
308  END DO
309  END IF
310  ELSE
311 *
312 * Multiply by inv(C**H).
313 *
314  IF ( capply ) THEN
315  DO i = 1, n
316  work( i ) = work( i ) * c( i )
317  END DO
318  END IF
319 *
320  IF ( notrans ) THEN
321  CALL zgbtrs( 'Conjugate transpose', n, kl, ku, 1, afb,
322  $ ldafb, ipiv, work, n, info )
323  ELSE
324  CALL zgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
325  $ ipiv, work, n, info )
326  END IF
327 *
328 * Multiply by R.
329 *
330  DO i = 1, n
331  work( i ) = work( i ) * rwork( i )
332  END DO
333  END IF
334  GO TO 10
335  END IF
336 *
337 * Compute the estimate of the reciprocal condition number.
338 *
339  IF( ainvnm .NE. 0.0d+0 )
340  $ zla_gbrcond_c = 1.0d+0 / ainvnm
341 *
342  RETURN
343 *
subroutine zgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
ZGBTRS
Definition: zgbtrs.f:140
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
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:135
double precision function zla_gbrcond_c(TRANS, N, KL, KU, AB, LDAB, AFB, LDAFB, IPIV, C, CAPPLY, INFO, WORK, RWORK)
ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general banded ma...
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: