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

◆ cgbt05()

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, A**T, or A**H, 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.

Definition at line 174 of file cgbt05.f.

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