LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zpbt05 ( character  UPLO,
integer  N,
integer  KD,
integer  NRHS,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
complex*16, dimension( ldb, * )  B,
integer  LDB,
complex*16, dimension( ldx, * )  X,
integer  LDX,
complex*16, dimension( ldxact, * )  XACT,
integer  LDXACT,
double precision, dimension( * )  FERR,
double precision, dimension( * )  BERR,
double precision, dimension( * )  RESLTS 
)

ZPBT05

Purpose:
 ZPBT05 tests the error bounds from iterative refinement for the
 computed solution to a system of equations A*X = B, where A is a
 Hermitian band matrix.

 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(A)*abs(X) +abs(b))_i )
             and NZ = max. number of nonzeros in any row of A, plus 1
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          Hermitian matrix A is stored.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[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]KD
          KD is INTEGER
          The number of super-diagonals of the matrix A if UPLO = 'U',
          or the number of sub-diagonals if UPLO = 'L'.  KD >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of columns of the matrices X, B, and XACT.
          NRHS >= 0.
[in]AB
          AB is COMPLEX*16 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]B
          B is COMPLEX*16 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*16 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*16 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 173 of file zpbt05.f.

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

Here is the caller graph for this function: