LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dpbt05()

subroutine dpbt05 ( character  uplo,
integer  n,
integer  kd,
integer  nrhs,
double precision, dimension( ldab, * )  ab,
integer  ldab,
double precision, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( ldx, * )  x,
integer  ldx,
double precision, dimension( ldxact, * )  xact,
integer  ldxact,
double precision, dimension( * )  ferr,
double precision, dimension( * )  berr,
double precision, dimension( * )  reslts 
)

DPBT05

Purpose:
 DPBT05 tests the error bounds from iterative refinement for the
 computed solution to a system of equations A*X = B, where A is a
 symmetric 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
          symmetric 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 DOUBLE PRECISION array, dimension (LDAB,N)
          The upper or lower triangle of the symmetric 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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.

Definition at line 169 of file dpbt05.f.

171*
172* -- LAPACK test routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 CHARACTER UPLO
178 INTEGER KD, LDAB, LDB, LDX, LDXACT, N, NRHS
179* ..
180* .. Array Arguments ..
181 DOUBLE PRECISION AB( LDAB, * ), B( LDB, * ), BERR( * ),
182 $ FERR( * ), RESLTS( * ), X( LDX, * ),
183 $ XACT( LDXACT, * )
184* ..
185*
186* =====================================================================
187*
188* .. Parameters ..
189 DOUBLE PRECISION ZERO, ONE
190 parameter( zero = 0.0d+0, one = 1.0d+0 )
191* ..
192* .. Local Scalars ..
193 LOGICAL UPPER
194 INTEGER I, IMAX, J, K, NZ
195 DOUBLE PRECISION AXBI, DIFF, EPS, ERRBND, OVFL, TMP, UNFL, XNORM
196* ..
197* .. External Functions ..
198 LOGICAL LSAME
199 INTEGER IDAMAX
200 DOUBLE PRECISION DLAMCH
201 EXTERNAL lsame, idamax, dlamch
202* ..
203* .. Intrinsic Functions ..
204 INTRINSIC abs, max, min
205* ..
206* .. Executable Statements ..
207*
208* Quick exit if N = 0 or NRHS = 0.
209*
210 IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
211 reslts( 1 ) = zero
212 reslts( 2 ) = zero
213 RETURN
214 END IF
215*
216 eps = dlamch( 'Epsilon' )
217 unfl = dlamch( 'Safe minimum' )
218 ovfl = one / unfl
219 upper = lsame( uplo, 'U' )
220 nz = 2*max( kd, n-1 ) + 1
221*
222* Test 1: Compute the maximum of
223* norm(X - XACT) / ( norm(X) * FERR )
224* over all the vectors X and XACT using the infinity-norm.
225*
226 errbnd = zero
227 DO 30 j = 1, nrhs
228 imax = idamax( n, x( 1, j ), 1 )
229 xnorm = max( abs( x( imax, j ) ), unfl )
230 diff = zero
231 DO 10 i = 1, n
232 diff = max( diff, abs( x( i, j )-xact( i, j ) ) )
233 10 CONTINUE
234*
235 IF( xnorm.GT.one ) THEN
236 GO TO 20
237 ELSE IF( diff.LE.ovfl*xnorm ) THEN
238 GO TO 20
239 ELSE
240 errbnd = one / eps
241 GO TO 30
242 END IF
243*
244 20 CONTINUE
245 IF( diff / xnorm.LE.ferr( j ) ) THEN
246 errbnd = max( errbnd, ( diff / xnorm ) / ferr( j ) )
247 ELSE
248 errbnd = one / eps
249 END IF
250 30 CONTINUE
251 reslts( 1 ) = errbnd
252*
253* Test 2: Compute the maximum of BERR / ( NZ*EPS + (*) ), where
254* (*) = NZ*UNFL / (min_i (abs(A)*abs(X) +abs(b))_i )
255*
256 DO 90 k = 1, nrhs
257 DO 80 i = 1, n
258 tmp = abs( b( i, k ) )
259 IF( upper ) THEN
260 DO 40 j = max( i-kd, 1 ), i
261 tmp = tmp + abs( ab( kd+1-i+j, i ) )*abs( x( j, k ) )
262 40 CONTINUE
263 DO 50 j = i + 1, min( i+kd, n )
264 tmp = tmp + abs( ab( kd+1+i-j, j ) )*abs( x( j, k ) )
265 50 CONTINUE
266 ELSE
267 DO 60 j = max( i-kd, 1 ), i - 1
268 tmp = tmp + abs( ab( 1+i-j, j ) )*abs( x( j, k ) )
269 60 CONTINUE
270 DO 70 j = i, min( i+kd, n )
271 tmp = tmp + abs( ab( 1+j-i, i ) )*abs( x( j, k ) )
272 70 CONTINUE
273 END IF
274 IF( i.EQ.1 ) THEN
275 axbi = tmp
276 ELSE
277 axbi = min( axbi, tmp )
278 END IF
279 80 CONTINUE
280 tmp = berr( k ) / ( nz*eps+nz*unfl / max( axbi, nz*unfl ) )
281 IF( k.EQ.1 ) THEN
282 reslts( 2 ) = tmp
283 ELSE
284 reslts( 2 ) = max( reslts( 2 ), tmp )
285 END IF
286 90 CONTINUE
287*
288 RETURN
289*
290* End of DPBT05
291*
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the caller graph for this function: