LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cgbrfs ( character  TRANS,
integer  N,
integer  KL,
integer  KU,
integer  NRHS,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldafb, * )  AFB,
integer  LDAFB,
integer, dimension( * )  IPIV,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldx, * )  X,
integer  LDX,
real, dimension( * )  FERR,
real, dimension( * )  BERR,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CGBRFS

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

Purpose:
 CGBRFS improves the computed solution to a system of linear
 equations when the coefficient matrix is banded, and provides
 error bounds and backward error estimates for the solution.
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)
[in]N
          N is INTEGER
          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]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrices B and X.  NRHS >= 0.
[in]AB
          AB is COMPLEX array, dimension (LDAB,N)
          The original band matrix A, stored 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 CGBTRF; for 1<=i<=N, row i of the
          matrix was interchanged with row IPIV(i).
[in]B
          B is COMPLEX array, dimension (LDB,NRHS)
          The right hand side matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in,out]X
          X is COMPLEX array, dimension (LDX,NRHS)
          On entry, the solution matrix X, as computed by CGBTRS.
          On exit, the improved solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]FERR
          FERR is REAL array, dimension (NRHS)
          The estimated forward error bound for each solution vector
          X(j) (the j-th column of the solution matrix X).
          If XTRUE is the true solution corresponding to X(j), FERR(j)
          is an estimated upper bound for the magnitude of the largest
          element in (X(j) - XTRUE) divided by the magnitude of the
          largest element in X(j).  The estimate is as reliable as
          the estimate for RCOND, and is almost always a slight
          overestimate of the true error.
[out]BERR
          BERR is REAL array, dimension (NRHS)
          The componentwise relative backward error of each solution
          vector X(j) (i.e., the smallest relative change in
          any element of A or B that makes X(j) an exact solution).
[out]WORK
          WORK is COMPLEX array, dimension (2*N)
[out]RWORK
          RWORK is REAL array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Internal Parameters:
  ITMAX is the maximum number of steps of iterative refinement.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 208 of file cgbrfs.f.

208 *
209 * -- LAPACK computational routine (version 3.4.0) --
210 * -- LAPACK is a software package provided by Univ. of Tennessee, --
211 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
212 * November 2011
213 *
214 * .. Scalar Arguments ..
215  CHARACTER trans
216  INTEGER info, kl, ku, ldab, ldafb, ldb, ldx, n, nrhs
217 * ..
218 * .. Array Arguments ..
219  INTEGER ipiv( * )
220  REAL berr( * ), ferr( * ), rwork( * )
221  COMPLEX ab( ldab, * ), afb( ldafb, * ), b( ldb, * ),
222  $ work( * ), x( ldx, * )
223 * ..
224 *
225 * =====================================================================
226 *
227 * .. Parameters ..
228  INTEGER itmax
229  parameter ( itmax = 5 )
230  REAL zero
231  parameter ( zero = 0.0e+0 )
232  COMPLEX cone
233  parameter ( cone = ( 1.0e+0, 0.0e+0 ) )
234  REAL two
235  parameter ( two = 2.0e+0 )
236  REAL three
237  parameter ( three = 3.0e+0 )
238 * ..
239 * .. Local Scalars ..
240  LOGICAL notran
241  CHARACTER transn, transt
242  INTEGER count, i, j, k, kase, kk, nz
243  REAL eps, lstres, s, safe1, safe2, safmin, xk
244  COMPLEX zdum
245 * ..
246 * .. Local Arrays ..
247  INTEGER isave( 3 )
248 * ..
249 * .. External Subroutines ..
250  EXTERNAL caxpy, ccopy, cgbmv, cgbtrs, clacn2, xerbla
251 * ..
252 * .. Intrinsic Functions ..
253  INTRINSIC abs, aimag, max, min, real
254 * ..
255 * .. External Functions ..
256  LOGICAL lsame
257  REAL slamch
258  EXTERNAL lsame, slamch
259 * ..
260 * .. Statement Functions ..
261  REAL cabs1
262 * ..
263 * .. Statement Function definitions ..
264  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
265 * ..
266 * .. Executable Statements ..
267 *
268 * Test the input parameters.
269 *
270  info = 0
271  notran = lsame( trans, 'N' )
272  IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
273  $ lsame( trans, 'C' ) ) THEN
274  info = -1
275  ELSE IF( n.LT.0 ) THEN
276  info = -2
277  ELSE IF( kl.LT.0 ) THEN
278  info = -3
279  ELSE IF( ku.LT.0 ) THEN
280  info = -4
281  ELSE IF( nrhs.LT.0 ) THEN
282  info = -5
283  ELSE IF( ldab.LT.kl+ku+1 ) THEN
284  info = -7
285  ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
286  info = -9
287  ELSE IF( ldb.LT.max( 1, n ) ) THEN
288  info = -12
289  ELSE IF( ldx.LT.max( 1, n ) ) THEN
290  info = -14
291  END IF
292  IF( info.NE.0 ) THEN
293  CALL xerbla( 'CGBRFS', -info )
294  RETURN
295  END IF
296 *
297 * Quick return if possible
298 *
299  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
300  DO 10 j = 1, nrhs
301  ferr( j ) = zero
302  berr( j ) = zero
303  10 CONTINUE
304  RETURN
305  END IF
306 *
307  IF( notran ) THEN
308  transn = 'N'
309  transt = 'C'
310  ELSE
311  transn = 'C'
312  transt = 'N'
313  END IF
314 *
315 * NZ = maximum number of nonzero elements in each row of A, plus 1
316 *
317  nz = min( kl+ku+2, n+1 )
318  eps = slamch( 'Epsilon' )
319  safmin = slamch( 'Safe minimum' )
320  safe1 = nz*safmin
321  safe2 = safe1 / eps
322 *
323 * Do for each right hand side
324 *
325  DO 140 j = 1, nrhs
326 *
327  count = 1
328  lstres = three
329  20 CONTINUE
330 *
331 * Loop until stopping criterion is satisfied.
332 *
333 * Compute residual R = B - op(A) * X,
334 * where op(A) = A, A**T, or A**H, depending on TRANS.
335 *
336  CALL ccopy( n, b( 1, j ), 1, work, 1 )
337  CALL cgbmv( trans, n, n, kl, ku, -cone, ab, ldab, x( 1, j ), 1,
338  $ cone, work, 1 )
339 *
340 * Compute componentwise relative backward error from formula
341 *
342 * max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
343 *
344 * where abs(Z) is the componentwise absolute value of the matrix
345 * or vector Z. If the i-th component of the denominator is less
346 * than SAFE2, then SAFE1 is added to the i-th components of the
347 * numerator and denominator before dividing.
348 *
349  DO 30 i = 1, n
350  rwork( i ) = cabs1( b( i, j ) )
351  30 CONTINUE
352 *
353 * Compute abs(op(A))*abs(X) + abs(B).
354 *
355  IF( notran ) THEN
356  DO 50 k = 1, n
357  kk = ku + 1 - k
358  xk = cabs1( x( k, j ) )
359  DO 40 i = max( 1, k-ku ), min( n, k+kl )
360  rwork( i ) = rwork( i ) + cabs1( ab( kk+i, k ) )*xk
361  40 CONTINUE
362  50 CONTINUE
363  ELSE
364  DO 70 k = 1, n
365  s = zero
366  kk = ku + 1 - k
367  DO 60 i = max( 1, k-ku ), min( n, k+kl )
368  s = s + cabs1( ab( kk+i, k ) )*cabs1( x( i, j ) )
369  60 CONTINUE
370  rwork( k ) = rwork( k ) + s
371  70 CONTINUE
372  END IF
373  s = zero
374  DO 80 i = 1, n
375  IF( rwork( i ).GT.safe2 ) THEN
376  s = max( s, cabs1( work( i ) ) / rwork( i ) )
377  ELSE
378  s = max( s, ( cabs1( work( i ) )+safe1 ) /
379  $ ( rwork( i )+safe1 ) )
380  END IF
381  80 CONTINUE
382  berr( j ) = s
383 *
384 * Test stopping criterion. Continue iterating if
385 * 1) The residual BERR(J) is larger than machine epsilon, and
386 * 2) BERR(J) decreased by at least a factor of 2 during the
387 * last iteration, and
388 * 3) At most ITMAX iterations tried.
389 *
390  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
391  $ count.LE.itmax ) THEN
392 *
393 * Update solution and try again.
394 *
395  CALL cgbtrs( trans, n, kl, ku, 1, afb, ldafb, ipiv, work, n,
396  $ info )
397  CALL caxpy( n, cone, work, 1, x( 1, j ), 1 )
398  lstres = berr( j )
399  count = count + 1
400  GO TO 20
401  END IF
402 *
403 * Bound error from formula
404 *
405 * norm(X - XTRUE) / norm(X) .le. FERR =
406 * norm( abs(inv(op(A)))*
407 * ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
408 *
409 * where
410 * norm(Z) is the magnitude of the largest component of Z
411 * inv(op(A)) is the inverse of op(A)
412 * abs(Z) is the componentwise absolute value of the matrix or
413 * vector Z
414 * NZ is the maximum number of nonzeros in any row of A, plus 1
415 * EPS is machine epsilon
416 *
417 * The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
418 * is incremented by SAFE1 if the i-th component of
419 * abs(op(A))*abs(X) + abs(B) is less than SAFE2.
420 *
421 * Use CLACN2 to estimate the infinity-norm of the matrix
422 * inv(op(A)) * diag(W),
423 * where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
424 *
425  DO 90 i = 1, n
426  IF( rwork( i ).GT.safe2 ) THEN
427  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i )
428  ELSE
429  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i ) +
430  $ safe1
431  END IF
432  90 CONTINUE
433 *
434  kase = 0
435  100 CONTINUE
436  CALL clacn2( n, work( n+1 ), work, ferr( j ), kase, isave )
437  IF( kase.NE.0 ) THEN
438  IF( kase.EQ.1 ) THEN
439 *
440 * Multiply by diag(W)*inv(op(A)**H).
441 *
442  CALL cgbtrs( transt, n, kl, ku, 1, afb, ldafb, ipiv,
443  $ work, n, info )
444  DO 110 i = 1, n
445  work( i ) = rwork( i )*work( i )
446  110 CONTINUE
447  ELSE
448 *
449 * Multiply by inv(op(A))*diag(W).
450 *
451  DO 120 i = 1, n
452  work( i ) = rwork( i )*work( i )
453  120 CONTINUE
454  CALL cgbtrs( transn, n, kl, ku, 1, afb, ldafb, ipiv,
455  $ work, n, info )
456  END IF
457  GO TO 100
458  END IF
459 *
460 * Normalize error.
461 *
462  lstres = zero
463  DO 130 i = 1, n
464  lstres = max( lstres, cabs1( x( i, j ) ) )
465  130 CONTINUE
466  IF( lstres.NE.zero )
467  $ ferr( j ) = ferr( j ) / lstres
468 *
469  140 CONTINUE
470 *
471  RETURN
472 *
473 * End of CGBRFS
474 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cgbmv(TRANS, M, N, KL, KU, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CGBMV
Definition: cgbmv.f:189
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53
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: