LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
double precision function zla_gbrcond_x ( 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,
complex*16, dimension( * )  X,
integer  INFO,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK 
)

ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrices.

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

Purpose:
    ZLA_GBRCOND_X Computes the infinity norm condition number of
    op(A) * diag(X) where X is a COMPLEX*16 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]X
          X is COMPLEX*16 array, dimension (N)
     The vector X in the formula op(A) * diag(X).
[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 157 of file zla_gbrcond_x.f.

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