LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cgbt05 ( character  TRANS,
integer  N,
integer  KL,
integer  KU,
integer  NRHS,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldx, * )  X,
integer  LDX,
complex, dimension( ldxact, * )  XACT,
integer  LDXACT,
real, dimension( * )  FERR,
real, dimension( * )  BERR,
real, dimension( * )  RESLTS 
)

CGBT05

Purpose:
 CGBT05 tests the error bounds from iterative refinement for the
 computed solution to a system of equations op(A)*X = B, where A is a
 general band matrix of order n with kl subdiagonals and ku
 superdiagonals and op(A) = A or A**T, depending on TRANS.

 RESLTS(1) = test of the error bound
           = norm(X - XACT) / ( norm(X) * FERR )

 A large value is returned if this ratio is not less than one.

 RESLTS(2) = residual from the iterative refinement routine
           = the maximum of BERR / ( NZ*EPS + (*) ), where
             (*) = NZ*UNFL / (min_i (abs(op(A))*abs(X) +abs(b))_i )
             and NZ = max. number of nonzeros in any row of A, plus 1
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 rows of the matrices X, B, and XACT, and 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 columns of the matrices X, B, and XACT.
          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]B
          B is COMPLEX array, dimension (LDB,NRHS)
          The right hand side vectors for the system of linear
          equations.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in]X
          X is COMPLEX array, dimension (LDX,NRHS)
          The computed solution vectors.  Each vector is stored as a
          column of the matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[in]XACT
          XACT is COMPLEX array, dimension (LDX,NRHS)
          The exact solution vectors.  Each vector is stored as a
          column of the matrix XACT.
[in]LDXACT
          LDXACT is INTEGER
          The leading dimension of the array XACT.  LDXACT >= max(1,N).
[in]FERR
          FERR is REAL array, dimension (NRHS)
          The estimated forward error bounds for each solution vector
          X.  If XTRUE is the true solution, FERR bounds the magnitude
          of the largest entry in (X - XTRUE) divided by the magnitude
          of the largest entry in X.
[in]BERR
          BERR is REAL array, dimension (NRHS)
          The componentwise relative backward error of each solution
          vector (i.e., the smallest relative change in any entry of A
          or B that makes X an exact solution).
[out]RESLTS
          RESLTS is REAL array, dimension (2)
          The maximum over the NRHS solution vectors of the ratios:
          RESLTS(1) = norm(X - XACT) / ( norm(X) * FERR )
          RESLTS(2) = BERR / ( NZ*EPS + (*) )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 178 of file cgbt05.f.

178 *
179 * -- LAPACK test routine (version 3.4.0) --
180 * -- LAPACK is a software package provided by Univ. of Tennessee, --
181 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
182 * November 2011
183 *
184 * .. Scalar Arguments ..
185  CHARACTER trans
186  INTEGER kl, ku, ldab, ldb, ldx, ldxact, n, nrhs
187 * ..
188 * .. Array Arguments ..
189  REAL berr( * ), ferr( * ), reslts( * )
190  COMPLEX ab( ldab, * ), b( ldb, * ), x( ldx, * ),
191  $ xact( ldxact, * )
192 * ..
193 *
194 * =====================================================================
195 *
196 * .. Parameters ..
197  REAL zero, one
198  parameter ( zero = 0.0e+0, one = 1.0e+0 )
199 * ..
200 * .. Local Scalars ..
201  LOGICAL notran
202  INTEGER i, imax, j, k, nz
203  REAL axbi, diff, eps, errbnd, ovfl, tmp, unfl, xnorm
204  COMPLEX zdum
205 * ..
206 * .. External Functions ..
207  LOGICAL lsame
208  INTEGER icamax
209  REAL slamch
210  EXTERNAL lsame, icamax, slamch
211 * ..
212 * .. Intrinsic Functions ..
213  INTRINSIC abs, aimag, max, min, real
214 * ..
215 * .. Statement Functions ..
216  REAL cabs1
217 * ..
218 * .. Statement Function definitions ..
219  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
220 * ..
221 * .. Executable Statements ..
222 *
223 * Quick exit if N = 0 or NRHS = 0.
224 *
225  IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
226  reslts( 1 ) = zero
227  reslts( 2 ) = zero
228  RETURN
229  END IF
230 *
231  eps = slamch( 'Epsilon' )
232  unfl = slamch( 'Safe minimum' )
233  ovfl = one / unfl
234  notran = lsame( trans, 'N' )
235  nz = min( kl+ku+2, n+1 )
236 *
237 * Test 1: Compute the maximum of
238 * norm(X - XACT) / ( norm(X) * FERR )
239 * over all the vectors X and XACT using the infinity-norm.
240 *
241  errbnd = zero
242  DO 30 j = 1, nrhs
243  imax = icamax( n, x( 1, j ), 1 )
244  xnorm = max( cabs1( x( imax, j ) ), unfl )
245  diff = zero
246  DO 10 i = 1, n
247  diff = max( diff, cabs1( x( i, j )-xact( i, j ) ) )
248  10 CONTINUE
249 *
250  IF( xnorm.GT.one ) THEN
251  GO TO 20
252  ELSE IF( diff.LE.ovfl*xnorm ) THEN
253  GO TO 20
254  ELSE
255  errbnd = one / eps
256  GO TO 30
257  END IF
258 *
259  20 CONTINUE
260  IF( diff / xnorm.LE.ferr( j ) ) THEN
261  errbnd = max( errbnd, ( diff / xnorm ) / ferr( j ) )
262  ELSE
263  errbnd = one / eps
264  END IF
265  30 CONTINUE
266  reslts( 1 ) = errbnd
267 *
268 * Test 2: Compute the maximum of BERR / ( NZ*EPS + (*) ), where
269 * (*) = NZ*UNFL / (min_i (abs(op(A))*abs(X) +abs(b))_i )
270 *
271  DO 70 k = 1, nrhs
272  DO 60 i = 1, n
273  tmp = cabs1( b( i, k ) )
274  IF( notran ) THEN
275  DO 40 j = max( i-kl, 1 ), min( i+ku, n )
276  tmp = tmp + cabs1( ab( ku+1+i-j, j ) )*
277  $ cabs1( x( j, k ) )
278  40 CONTINUE
279  ELSE
280  DO 50 j = max( i-ku, 1 ), min( i+kl, n )
281  tmp = tmp + cabs1( ab( ku+1+j-i, i ) )*
282  $ cabs1( x( j, k ) )
283  50 CONTINUE
284  END IF
285  IF( i.EQ.1 ) THEN
286  axbi = tmp
287  ELSE
288  axbi = min( axbi, tmp )
289  END IF
290  60 CONTINUE
291  tmp = berr( k ) / ( nz*eps+nz*unfl / max( axbi, nz*unfl ) )
292  IF( k.EQ.1 ) THEN
293  reslts( 2 ) = tmp
294  ELSE
295  reslts( 2 ) = max( reslts( 2 ), tmp )
296  END IF
297  70 CONTINUE
298 *
299  RETURN
300 *
301 * End of CGBT05
302 *
integer function icamax(N, CX, INCX)
ICAMAX
Definition: icamax.f:53
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the caller graph for this function: