LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
real function sla_gbrcond ( character  TRANS,
integer  N,
integer  KL,
integer  KU,
real, dimension( ldab, * )  AB,
integer  LDAB,
real, dimension( ldafb, * )  AFB,
integer  LDAFB,
integer, dimension( * )  IPIV,
integer  CMODE,
real, dimension( * )  C,
integer  INFO,
real, dimension( * )  WORK,
integer, dimension( * )  IWORK 
)

SLA_GBRCOND estimates the Skeel condition number for a general banded matrix.

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

Purpose:
    SLA_GBRCOND Estimates the Skeel condition number of  op(A) * op2(C)
    where op2 is determined by CMODE as follows
    CMODE =  1    op2(C) = C
    CMODE =  0    op2(C) = I
    CMODE = -1    op2(C) = inv(C)
    The Skeel condition number  cond(A) = norminf( |inv(A)||A| )
    is computed by computing scaling factors R such that
    diag(R)*A*op2(C) is row equilibrated and computing the standard
    infinity-norm condition number.
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 REAL 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 REAL array, dimension (LDAFB,N)
     Details of the LU factorization of the band matrix A, as
     computed by SGBTRF.  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 SGBTRF; row i of the matrix was interchanged
     with row IPIV(i).
[in]CMODE
          CMODE is INTEGER
     Determines op2(C) in the formula op(A) * op2(C) as follows:
     CMODE =  1    op2(C) = C
     CMODE =  0    op2(C) = I
     CMODE = -1    op2(C) = inv(C)
[in]C
          C is REAL array, dimension (N)
     The vector C in the formula op(A) * op2(C).
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[in]WORK
          WORK is REAL array, dimension (5*N).
     Workspace.
[in]IWORK
          IWORK is INTEGER array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 170 of file sla_gbrcond.f.

170 *
171 * -- LAPACK computational routine (version 3.4.2) --
172 * -- LAPACK is a software package provided by Univ. of Tennessee, --
173 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
174 * September 2012
175 *
176 * .. Scalar Arguments ..
177  CHARACTER trans
178  INTEGER n, ldab, ldafb, info, kl, ku, cmode
179 * ..
180 * .. Array Arguments ..
181  INTEGER iwork( * ), ipiv( * )
182  REAL ab( ldab, * ), afb( ldafb, * ), work( * ),
183  $ c( * )
184 * ..
185 *
186 * =====================================================================
187 *
188 * .. Local Scalars ..
189  LOGICAL notrans
190  INTEGER kase, i, j, kd, ke
191  REAL ainvnm, tmp
192 * ..
193 * .. Local Arrays ..
194  INTEGER isave( 3 )
195 * ..
196 * .. External Functions ..
197  LOGICAL lsame
198  EXTERNAL lsame
199 * ..
200 * .. External Subroutines ..
201  EXTERNAL slacn2, sgbtrs, xerbla
202 * ..
203 * .. Intrinsic Functions ..
204  INTRINSIC abs, max
205 * ..
206 * .. Executable Statements ..
207 *
208  sla_gbrcond = 0.0
209 *
210  info = 0
211  notrans = lsame( trans, 'N' )
212  IF ( .NOT. notrans .AND. .NOT. lsame(trans, 'T')
213  $ .AND. .NOT. 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( 'SLA_GBRCOND', -info )
228  RETURN
229  END IF
230  IF( n.EQ.0 ) THEN
231  sla_gbrcond = 1.0
232  RETURN
233  END IF
234 *
235 * Compute the equilibration matrix R such that
236 * inv(R)*A*C has unit 1-norm.
237 *
238  kd = ku + 1
239  ke = kl + 1
240  IF ( notrans ) THEN
241  DO i = 1, n
242  tmp = 0.0
243  IF ( cmode .EQ. 1 ) THEN
244  DO j = max( i-kl, 1 ), min( i+ku, n )
245  tmp = tmp + abs( ab( kd+i-j, j ) * c( j ) )
246  END DO
247  ELSE IF ( cmode .EQ. 0 ) THEN
248  DO j = max( i-kl, 1 ), min( i+ku, n )
249  tmp = tmp + abs( ab( kd+i-j, j ) )
250  END DO
251  ELSE
252  DO j = max( i-kl, 1 ), min( i+ku, n )
253  tmp = tmp + abs( ab( kd+i-j, j ) / c( j ) )
254  END DO
255  END IF
256  work( 2*n+i ) = tmp
257  END DO
258  ELSE
259  DO i = 1, n
260  tmp = 0.0
261  IF ( cmode .EQ. 1 ) THEN
262  DO j = max( i-kl, 1 ), min( i+ku, n )
263  tmp = tmp + abs( ab( ke-i+j, i ) * c( j ) )
264  END DO
265  ELSE IF ( cmode .EQ. 0 ) THEN
266  DO j = max( i-kl, 1 ), min( i+ku, n )
267  tmp = tmp + abs( ab( ke-i+j, i ) )
268  END DO
269  ELSE
270  DO j = max( i-kl, 1 ), min( i+ku, n )
271  tmp = tmp + abs( ab( ke-i+j, i ) / c( j ) )
272  END DO
273  END IF
274  work( 2*n+i ) = tmp
275  END DO
276  END IF
277 *
278 * Estimate the norm of inv(op(A)).
279 *
280  ainvnm = 0.0
281 
282  kase = 0
283  10 CONTINUE
284  CALL slacn2( n, work( n+1 ), work, iwork, ainvnm, kase, isave )
285  IF( kase.NE.0 ) THEN
286  IF( kase.EQ.2 ) THEN
287 *
288 * Multiply by R.
289 *
290  DO i = 1, n
291  work( i ) = work( i ) * work( 2*n+i )
292  END DO
293 
294  IF ( notrans ) THEN
295  CALL sgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
296  $ ipiv, work, n, info )
297  ELSE
298  CALL sgbtrs( 'Transpose', n, kl, ku, 1, afb, ldafb, ipiv,
299  $ work, n, info )
300  END IF
301 *
302 * Multiply by inv(C).
303 *
304  IF ( cmode .EQ. 1 ) THEN
305  DO i = 1, n
306  work( i ) = work( i ) / c( i )
307  END DO
308  ELSE IF ( cmode .EQ. -1 ) THEN
309  DO i = 1, n
310  work( i ) = work( i ) * c( i )
311  END DO
312  END IF
313  ELSE
314 *
315 * Multiply by inv(C**T).
316 *
317  IF ( cmode .EQ. 1 ) THEN
318  DO i = 1, n
319  work( i ) = work( i ) / c( i )
320  END DO
321  ELSE IF ( cmode .EQ. -1 ) THEN
322  DO i = 1, n
323  work( i ) = work( i ) * c( i )
324  END DO
325  END IF
326 
327  IF ( notrans ) THEN
328  CALL sgbtrs( 'Transpose', n, kl, ku, 1, afb, ldafb, ipiv,
329  $ work, n, info )
330  ELSE
331  CALL sgbtrs( 'No transpose', n, kl, ku, 1, afb, ldafb,
332  $ ipiv, work, n, info )
333  END IF
334 *
335 * Multiply by R.
336 *
337  DO i = 1, n
338  work( i ) = work( i ) * work( 2*n+i )
339  END DO
340  END IF
341  GO TO 10
342  END IF
343 *
344 * Compute the estimate of the reciprocal condition number.
345 *
346  IF( ainvnm .NE. 0.0 )
347  $ sla_gbrcond = ( 1.0 / ainvnm )
348 *
349  RETURN
350 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function sla_gbrcond(TRANS, N, KL, KU, AB, LDAB, AFB, LDAFB, IPIV, CMODE, C, INFO, WORK, IWORK)
SLA_GBRCOND estimates the Skeel condition number for a general banded matrix.
Definition: sla_gbrcond.f:170
subroutine slacn2(N, V, X, ISGN, EST, KASE, ISAVE)
SLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: slacn2.f:138
subroutine sgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
SGBTRS
Definition: sgbtrs.f:140
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: