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

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

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

Purpose:
    CLA_GBRCOND_C Computes the infinity norm condition number of
    op(A) * inv(diag(C)) where C is a REAL 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]C
          C is REAL 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 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 163 of file cla_gbrcond_c.f.

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