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

SGBT05

Purpose:
 SGBT05 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 REAL 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 REAL 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 REAL 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 REAL 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 sgbt05.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 ab( ldab, * ), b( ldb, * ), berr( * ),
190  $ ferr( * ), reslts( * ), 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 * ..
205 * .. External Functions ..
206  LOGICAL lsame
207  INTEGER isamax
208  REAL slamch
209  EXTERNAL lsame, isamax, slamch
210 * ..
211 * .. Intrinsic Functions ..
212  INTRINSIC abs, max, min
213 * ..
214 * .. Executable Statements ..
215 *
216 * Quick exit if N = 0 or NRHS = 0.
217 *
218  IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
219  reslts( 1 ) = zero
220  reslts( 2 ) = zero
221  RETURN
222  END IF
223 *
224  eps = slamch( 'Epsilon' )
225  unfl = slamch( 'Safe minimum' )
226  ovfl = one / unfl
227  notran = lsame( trans, 'N' )
228  nz = min( kl+ku+2, n+1 )
229 *
230 * Test 1: Compute the maximum of
231 * norm(X - XACT) / ( norm(X) * FERR )
232 * over all the vectors X and XACT using the infinity-norm.
233 *
234  errbnd = zero
235  DO 30 j = 1, nrhs
236  imax = isamax( n, x( 1, j ), 1 )
237  xnorm = max( abs( x( imax, j ) ), unfl )
238  diff = zero
239  DO 10 i = 1, n
240  diff = max( diff, abs( x( i, j )-xact( i, j ) ) )
241  10 CONTINUE
242 *
243  IF( xnorm.GT.one ) THEN
244  GO TO 20
245  ELSE IF( diff.LE.ovfl*xnorm ) THEN
246  GO TO 20
247  ELSE
248  errbnd = one / eps
249  GO TO 30
250  END IF
251 *
252  20 CONTINUE
253  IF( diff / xnorm.LE.ferr( j ) ) THEN
254  errbnd = max( errbnd, ( diff / xnorm ) / ferr( j ) )
255  ELSE
256  errbnd = one / eps
257  END IF
258  30 CONTINUE
259  reslts( 1 ) = errbnd
260 *
261 * Test 2: Compute the maximum of BERR / ( NZ*EPS + (*) ), where
262 * (*) = NZ*UNFL / (min_i (abs(op(A))*abs(X) +abs(b))_i )
263 *
264  DO 70 k = 1, nrhs
265  DO 60 i = 1, n
266  tmp = abs( b( i, k ) )
267  IF( notran ) THEN
268  DO 40 j = max( i-kl, 1 ), min( i+ku, n )
269  tmp = tmp + abs( ab( ku+1+i-j, j ) )*abs( x( j, k ) )
270  40 CONTINUE
271  ELSE
272  DO 50 j = max( i-ku, 1 ), min( i+kl, n )
273  tmp = tmp + abs( ab( ku+1+j-i, i ) )*abs( x( j, k ) )
274  50 CONTINUE
275  END IF
276  IF( i.EQ.1 ) THEN
277  axbi = tmp
278  ELSE
279  axbi = min( axbi, tmp )
280  END IF
281  60 CONTINUE
282  tmp = berr( k ) / ( nz*eps+nz*unfl / max( axbi, nz*unfl ) )
283  IF( k.EQ.1 ) THEN
284  reslts( 2 ) = tmp
285  ELSE
286  reslts( 2 ) = max( reslts( 2 ), tmp )
287  END IF
288  70 CONTINUE
289 *
290  RETURN
291 *
292 * End of SGBT05
293 *
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.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: