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

◆ zppt05()

subroutine zppt05 ( character  uplo,
integer  n,
integer  nrhs,
complex*16, dimension( * )  ap,
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 
)

ZPPT05

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

 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 / ( (n+1)*EPS + (*) ), where
             (*) = (n+1)*UNFL / (min_i (abs(A)*abs(X) +abs(b))_i )
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]NRHS
          NRHS is INTEGER
          The number of columns of the matrices X, B, and XACT.
          NRHS >= 0.
[in]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          The upper or lower triangle of the Hermitian matrix A, packed
          columnwise in a linear array.  The j-th column of A is stored
          in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
[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 / ( (n+1)*EPS + (*) )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 155 of file zppt05.f.

157*
158* -- LAPACK test routine --
159* -- LAPACK is a software package provided by Univ. of Tennessee, --
160* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
161*
162* .. Scalar Arguments ..
163 CHARACTER UPLO
164 INTEGER LDB, LDX, LDXACT, N, NRHS
165* ..
166* .. Array Arguments ..
167 DOUBLE PRECISION BERR( * ), FERR( * ), RESLTS( * )
168 COMPLEX*16 AP( * ), B( LDB, * ), X( LDX, * ),
169 $ XACT( LDXACT, * )
170* ..
171*
172* =====================================================================
173*
174* .. Parameters ..
175 DOUBLE PRECISION ZERO, ONE
176 parameter( zero = 0.0d+0, one = 1.0d+0 )
177* ..
178* .. Local Scalars ..
179 LOGICAL UPPER
180 INTEGER I, IMAX, J, JC, K
181 DOUBLE PRECISION AXBI, DIFF, EPS, ERRBND, OVFL, TMP, UNFL, XNORM
182 COMPLEX*16 ZDUM
183* ..
184* .. External Functions ..
185 LOGICAL LSAME
186 INTEGER IZAMAX
187 DOUBLE PRECISION DLAMCH
188 EXTERNAL lsame, izamax, dlamch
189* ..
190* .. Intrinsic Functions ..
191 INTRINSIC abs, dble, dimag, max, min
192* ..
193* .. Statement Functions ..
194 DOUBLE PRECISION CABS1
195* ..
196* .. Statement Function definitions ..
197 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
198* ..
199* .. Executable Statements ..
200*
201* Quick exit if N = 0 or NRHS = 0.
202*
203 IF( n.LE.0 .OR. nrhs.LE.0 ) THEN
204 reslts( 1 ) = zero
205 reslts( 2 ) = zero
206 RETURN
207 END IF
208*
209 eps = dlamch( 'Epsilon' )
210 unfl = dlamch( 'Safe minimum' )
211 ovfl = one / unfl
212 upper = lsame( uplo, 'U' )
213*
214* Test 1: Compute the maximum of
215* norm(X - XACT) / ( norm(X) * FERR )
216* over all the vectors X and XACT using the infinity-norm.
217*
218 errbnd = zero
219 DO 30 j = 1, nrhs
220 imax = izamax( n, x( 1, j ), 1 )
221 xnorm = max( cabs1( x( imax, j ) ), unfl )
222 diff = zero
223 DO 10 i = 1, n
224 diff = max( diff, cabs1( x( i, j )-xact( i, j ) ) )
225 10 CONTINUE
226*
227 IF( xnorm.GT.one ) THEN
228 GO TO 20
229 ELSE IF( diff.LE.ovfl*xnorm ) THEN
230 GO TO 20
231 ELSE
232 errbnd = one / eps
233 GO TO 30
234 END IF
235*
236 20 CONTINUE
237 IF( diff / xnorm.LE.ferr( j ) ) THEN
238 errbnd = max( errbnd, ( diff / xnorm ) / ferr( j ) )
239 ELSE
240 errbnd = one / eps
241 END IF
242 30 CONTINUE
243 reslts( 1 ) = errbnd
244*
245* Test 2: Compute the maximum of BERR / ( (n+1)*EPS + (*) ), where
246* (*) = (n+1)*UNFL / (min_i (abs(A)*abs(X) +abs(b))_i )
247*
248 DO 90 k = 1, nrhs
249 DO 80 i = 1, n
250 tmp = cabs1( b( i, k ) )
251 IF( upper ) THEN
252 jc = ( ( i-1 )*i ) / 2
253 DO 40 j = 1, i - 1
254 tmp = tmp + cabs1( ap( jc+j ) )*cabs1( x( j, k ) )
255 40 CONTINUE
256 tmp = tmp + abs( dble( ap( jc+i ) ) )*cabs1( x( i, k ) )
257 jc = jc + i + i
258 DO 50 j = i + 1, n
259 tmp = tmp + cabs1( ap( jc ) )*cabs1( x( j, k ) )
260 jc = jc + j
261 50 CONTINUE
262 ELSE
263 jc = i
264 DO 60 j = 1, i - 1
265 tmp = tmp + cabs1( ap( jc ) )*cabs1( x( j, k ) )
266 jc = jc + n - j
267 60 CONTINUE
268 tmp = tmp + abs( dble( ap( jc ) ) )*cabs1( x( i, k ) )
269 DO 70 j = i + 1, n
270 tmp = tmp + cabs1( ap( jc+j-i ) )*cabs1( x( j, k ) )
271 70 CONTINUE
272 END IF
273 IF( i.EQ.1 ) THEN
274 axbi = tmp
275 ELSE
276 axbi = min( axbi, tmp )
277 END IF
278 80 CONTINUE
279 tmp = berr( k ) / ( ( n+1 )*eps+( n+1 )*unfl /
280 $ max( axbi, ( n+1 )*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 ZPPT05
291*
integer function izamax(n, zx, incx)
IZAMAX
Definition izamax.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: