LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cpbrfs ( character  UPLO,
integer  N,
integer  KD,
integer  NRHS,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldafb, * )  AFB,
integer  LDAFB,
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 
)

CPBRFS

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

Purpose:
 CPBRFS improves the computed solution to a system of linear
 equations when the coefficient matrix is Hermitian positive definite
 and banded, and provides error bounds and backward error estimates
 for the solution.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 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 upper or lower triangle of the Hermitian band matrix A,
          stored in the first KD+1 rows of the array.  The j-th column
          of A is stored in the j-th column of the array AB as follows:
          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD+1.
[in]AFB
          AFB is COMPLEX array, dimension (LDAFB,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**H*U or A = L*L**H of the band matrix A as computed by
          CPBTRF, in the same storage format as A (see AB).
[in]LDAFB
          LDAFB is INTEGER
          The leading dimension of the array AFB.  LDAFB >= KD+1.
[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 CPBTRS.
          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 191 of file cpbrfs.f.

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