LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
real function cla_gbrcond_x ( character  TRANS,
integer  N,
integer  KL,
integer  KU,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldafb, * )  AFB,
integer  LDAFB,
integer, dimension( * )  IPIV,
complex, dimension( * )  X,
integer  INFO,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK 
)

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

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

Purpose:
    CLA_GBRCOND_X Computes the infinity norm condition number of
    op(A) * diag(X) where X is a COMPLEX 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 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 array, dimension (LDAFB,N)
     Details of the LU factorization of the band matrix A, as
     computed by CGBTRF.  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 CGBTRF; row i of the matrix was interchanged
     with row IPIV(i).
[in]X
          X is COMPLEX 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 array, dimension (2*N).
     Workspace.
[in]RWORK
          RWORK is REAL array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 155 of file cla_gbrcond_x.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: