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

◆ zchkpt()

subroutine zchkpt ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
logical  tsterr,
complex*16, dimension( * )  a,
double precision, dimension( * )  d,
complex*16, dimension( * )  e,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  xact,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer  nout 
)

ZCHKPT

Purpose:
 ZCHKPT tests ZPTTRF, -TRS, -RFS, and -CON
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*2)
[out]D
          D is DOUBLE PRECISION array, dimension (NMAX*2)
[out]E
          E is COMPLEX*16 array, dimension (NMAX*2)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NSMAX))
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 145 of file zchkpt.f.

147*
148* -- LAPACK test routine --
149* -- LAPACK is a software package provided by Univ. of Tennessee, --
150* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151*
152* .. Scalar Arguments ..
153 LOGICAL TSTERR
154 INTEGER NN, NNS, NOUT
155 DOUBLE PRECISION THRESH
156* ..
157* .. Array Arguments ..
158 LOGICAL DOTYPE( * )
159 INTEGER NSVAL( * ), NVAL( * )
160 DOUBLE PRECISION D( * ), RWORK( * )
161 COMPLEX*16 A( * ), B( * ), E( * ), WORK( * ), X( * ),
162 $ XACT( * )
163* ..
164*
165* =====================================================================
166*
167* .. Parameters ..
168 DOUBLE PRECISION ONE, ZERO
169 parameter( one = 1.0d+0, zero = 0.0d+0 )
170 INTEGER NTYPES
171 parameter( ntypes = 12 )
172 INTEGER NTESTS
173 parameter( ntests = 7 )
174* ..
175* .. Local Scalars ..
176 LOGICAL ZEROT
177 CHARACTER DIST, TYPE, UPLO
178 CHARACTER*3 PATH
179 INTEGER I, IA, IMAT, IN, INFO, IRHS, IUPLO, IX, IZERO,
180 $ J, K, KL, KU, LDA, MODE, N, NERRS, NFAIL,
181 $ NIMAT, NRHS, NRUN
182 DOUBLE PRECISION AINVNM, ANORM, COND, DMAX, RCOND, RCONDC
183* ..
184* .. Local Arrays ..
185 CHARACTER UPLOS( 2 )
186 INTEGER ISEED( 4 ), ISEEDY( 4 )
187 DOUBLE PRECISION RESULT( NTESTS )
188 COMPLEX*16 Z( 3 )
189* ..
190* .. External Functions ..
191 INTEGER IDAMAX
192 DOUBLE PRECISION DGET06, DZASUM, ZLANHT
193 EXTERNAL idamax, dget06, dzasum, zlanht
194* ..
195* .. External Subroutines ..
196 EXTERNAL alaerh, alahd, alasum, dcopy, dlarnv, dscal,
200* ..
201* .. Intrinsic Functions ..
202 INTRINSIC abs, dble, max
203* ..
204* .. Scalars in Common ..
205 LOGICAL LERR, OK
206 CHARACTER*32 SRNAMT
207 INTEGER INFOT, NUNIT
208* ..
209* .. Common blocks ..
210 COMMON / infoc / infot, nunit, ok, lerr
211 COMMON / srnamc / srnamt
212* ..
213* .. Data statements ..
214 DATA iseedy / 0, 0, 0, 1 / , uplos / 'U', 'L' /
215* ..
216* .. Executable Statements ..
217*
218 path( 1: 1 ) = 'Zomplex precision'
219 path( 2: 3 ) = 'PT'
220 nrun = 0
221 nfail = 0
222 nerrs = 0
223 DO 10 i = 1, 4
224 iseed( i ) = iseedy( i )
225 10 CONTINUE
226*
227* Test the error exits
228*
229 IF( tsterr )
230 $ CALL zerrgt( path, nout )
231 infot = 0
232*
233 DO 120 in = 1, nn
234*
235* Do for each value of N in NVAL.
236*
237 n = nval( in )
238 lda = max( 1, n )
239 nimat = ntypes
240 IF( n.LE.0 )
241 $ nimat = 1
242*
243 DO 110 imat = 1, nimat
244*
245* Do the tests only if DOTYPE( IMAT ) is true.
246*
247 IF( n.GT.0 .AND. .NOT.dotype( imat ) )
248 $ GO TO 110
249*
250* Set up parameters with ZLATB4.
251*
252 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
253 $ COND, DIST )
254*
255 zerot = imat.GE.8 .AND. imat.LE.10
256 IF( imat.LE.6 ) THEN
257*
258* Type 1-6: generate a Hermitian tridiagonal matrix of
259* known condition number in lower triangular band storage.
260*
261 srnamt = 'ZLATMS'
262 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE, COND,
263 $ ANORM, KL, KU, 'B', A, 2, WORK, INFO )
264*
265* Check the error code from ZLATMS.
266*
267 IF( info.NE.0 ) THEN
268 CALL alaerh( path, 'ZLATMS', info, 0, ' ', n, n, kl,
269 $ ku, -1, imat, nfail, nerrs, nout )
270 GO TO 110
271 END IF
272 izero = 0
273*
274* Copy the matrix to D and E.
275*
276 ia = 1
277 DO 20 i = 1, n - 1
278 d( i ) = dble( a( ia ) )
279 e( i ) = a( ia+1 )
280 ia = ia + 2
281 20 CONTINUE
282 IF( n.GT.0 )
283 $ d( n ) = dble( a( ia ) )
284 ELSE
285*
286* Type 7-12: generate a diagonally dominant matrix with
287* unknown condition number in the vectors D and E.
288*
289 IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
290*
291* Let E be complex, D real, with values from [-1,1].
292*
293 CALL dlarnv( 2, iseed, n, d )
294 CALL zlarnv( 2, iseed, n-1, e )
295*
296* Make the tridiagonal matrix diagonally dominant.
297*
298 IF( n.EQ.1 ) THEN
299 d( 1 ) = abs( d( 1 ) )
300 ELSE
301 d( 1 ) = abs( d( 1 ) ) + abs( e( 1 ) )
302 d( n ) = abs( d( n ) ) + abs( e( n-1 ) )
303 DO 30 i = 2, n - 1
304 d( i ) = abs( d( i ) ) + abs( e( i ) ) +
305 $ abs( e( i-1 ) )
306 30 CONTINUE
307 END IF
308*
309* Scale D and E so the maximum element is ANORM.
310*
311 ix = idamax( n, d, 1 )
312 dmax = d( ix )
313 CALL dscal( n, anorm / dmax, d, 1 )
314 CALL zdscal( n-1, anorm / dmax, e, 1 )
315*
316 ELSE IF( izero.GT.0 ) THEN
317*
318* Reuse the last matrix by copying back the zeroed out
319* elements.
320*
321 IF( izero.EQ.1 ) THEN
322 d( 1 ) = dble( z( 2 ) )
323 IF( n.GT.1 )
324 $ e( 1 ) = z( 3 )
325 ELSE IF( izero.EQ.n ) THEN
326 e( n-1 ) = z( 1 )
327 d( n ) = dble( z( 2 ) )
328 ELSE
329 e( izero-1 ) = z( 1 )
330 d( izero ) = dble( z( 2 ) )
331 e( izero ) = z( 3 )
332 END IF
333 END IF
334*
335* For types 8-10, set one row and column of the matrix to
336* zero.
337*
338 izero = 0
339 IF( imat.EQ.8 ) THEN
340 izero = 1
341 z( 2 ) = d( 1 )
342 d( 1 ) = zero
343 IF( n.GT.1 ) THEN
344 z( 3 ) = e( 1 )
345 e( 1 ) = zero
346 END IF
347 ELSE IF( imat.EQ.9 ) THEN
348 izero = n
349 IF( n.GT.1 ) THEN
350 z( 1 ) = e( n-1 )
351 e( n-1 ) = zero
352 END IF
353 z( 2 ) = d( n )
354 d( n ) = zero
355 ELSE IF( imat.EQ.10 ) THEN
356 izero = ( n+1 ) / 2
357 IF( izero.GT.1 ) THEN
358 z( 1 ) = e( izero-1 )
359 z( 3 ) = e( izero )
360 e( izero-1 ) = zero
361 e( izero ) = zero
362 END IF
363 z( 2 ) = d( izero )
364 d( izero ) = zero
365 END IF
366 END IF
367*
368 CALL dcopy( n, d, 1, d( n+1 ), 1 )
369 IF( n.GT.1 )
370 $ CALL zcopy( n-1, e, 1, e( n+1 ), 1 )
371*
372*+ TEST 1
373* Factor A as L*D*L' and compute the ratio
374* norm(L*D*L' - A) / (n * norm(A) * EPS )
375*
376 CALL zpttrf( n, d( n+1 ), e( n+1 ), info )
377*
378* Check error code from ZPTTRF.
379*
380 IF( info.NE.izero ) THEN
381 CALL alaerh( path, 'ZPTTRF', info, izero, ' ', n, n, -1,
382 $ -1, -1, imat, nfail, nerrs, nout )
383 GO TO 110
384 END IF
385*
386 IF( info.GT.0 ) THEN
387 rcondc = zero
388 GO TO 100
389 END IF
390*
391 CALL zptt01( n, d, e, d( n+1 ), e( n+1 ), work,
392 $ result( 1 ) )
393*
394* Print the test ratio if greater than or equal to THRESH.
395*
396 IF( result( 1 ).GE.thresh ) THEN
397 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
398 $ CALL alahd( nout, path )
399 WRITE( nout, fmt = 9999 )n, imat, 1, result( 1 )
400 nfail = nfail + 1
401 END IF
402 nrun = nrun + 1
403*
404* Compute RCONDC = 1 / (norm(A) * norm(inv(A))
405*
406* Compute norm(A).
407*
408 anorm = zlanht( '1', n, d, e )
409*
410* Use ZPTTRS to solve for one column at a time of inv(A),
411* computing the maximum column sum as we go.
412*
413 ainvnm = zero
414 DO 50 i = 1, n
415 DO 40 j = 1, n
416 x( j ) = zero
417 40 CONTINUE
418 x( i ) = one
419 CALL zpttrs( 'Lower', n, 1, d( n+1 ), e( n+1 ), x, lda,
420 $ info )
421 ainvnm = max( ainvnm, dzasum( n, x, 1 ) )
422 50 CONTINUE
423 rcondc = one / max( one, anorm*ainvnm )
424*
425 DO 90 irhs = 1, nns
426 nrhs = nsval( irhs )
427*
428* Generate NRHS random solution vectors.
429*
430 ix = 1
431 DO 60 j = 1, nrhs
432 CALL zlarnv( 2, iseed, n, xact( ix ) )
433 ix = ix + lda
434 60 CONTINUE
435*
436 DO 80 iuplo = 1, 2
437*
438* Do first for UPLO = 'U', then for UPLO = 'L'.
439*
440 uplo = uplos( iuplo )
441*
442* Set the right hand side.
443*
444 CALL zlaptm( uplo, n, nrhs, one, d, e, xact, lda,
445 $ zero, b, lda )
446*
447*+ TEST 2
448* Solve A*x = b and compute the residual.
449*
450 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
451 CALL zpttrs( uplo, n, nrhs, d( n+1 ), e( n+1 ), x,
452 $ lda, info )
453*
454* Check error code from ZPTTRS.
455*
456 IF( info.NE.0 )
457 $ CALL alaerh( path, 'ZPTTRS', info, 0, uplo, n, n,
458 $ -1, -1, nrhs, imat, nfail, nerrs,
459 $ nout )
460*
461 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
462 CALL zptt02( uplo, n, nrhs, d, e, x, lda, work, lda,
463 $ result( 2 ) )
464*
465*+ TEST 3
466* Check solution from generated exact solution.
467*
468 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
469 $ result( 3 ) )
470*
471*+ TESTS 4, 5, and 6
472* Use iterative refinement to improve the solution.
473*
474 srnamt = 'ZPTRFS'
475 CALL zptrfs( uplo, n, nrhs, d, e, d( n+1 ), e( n+1 ),
476 $ b, lda, x, lda, rwork, rwork( nrhs+1 ),
477 $ work, rwork( 2*nrhs+1 ), info )
478*
479* Check error code from ZPTRFS.
480*
481 IF( info.NE.0 )
482 $ CALL alaerh( path, 'ZPTRFS', info, 0, uplo, n, n,
483 $ -1, -1, nrhs, imat, nfail, nerrs,
484 $ nout )
485*
486 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
487 $ result( 4 ) )
488 CALL zptt05( n, nrhs, d, e, b, lda, x, lda, xact, lda,
489 $ rwork, rwork( nrhs+1 ), result( 5 ) )
490*
491* Print information about the tests that did not pass the
492* threshold.
493*
494 DO 70 k = 2, 6
495 IF( result( k ).GE.thresh ) THEN
496 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
497 $ CALL alahd( nout, path )
498 WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat,
499 $ k, result( k )
500 nfail = nfail + 1
501 END IF
502 70 CONTINUE
503 nrun = nrun + 5
504*
505 80 CONTINUE
506 90 CONTINUE
507*
508*+ TEST 7
509* Estimate the reciprocal of the condition number of the
510* matrix.
511*
512 100 CONTINUE
513 srnamt = 'ZPTCON'
514 CALL zptcon( n, d( n+1 ), e( n+1 ), anorm, rcond, rwork,
515 $ info )
516*
517* Check error code from ZPTCON.
518*
519 IF( info.NE.0 )
520 $ CALL alaerh( path, 'ZPTCON', info, 0, ' ', n, n, -1, -1,
521 $ -1, imat, nfail, nerrs, nout )
522*
523 result( 7 ) = dget06( rcond, rcondc )
524*
525* Print the test ratio if greater than or equal to THRESH.
526*
527 IF( result( 7 ).GE.thresh ) THEN
528 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
529 $ CALL alahd( nout, path )
530 WRITE( nout, fmt = 9999 )n, imat, 7, result( 7 )
531 nfail = nfail + 1
532 END IF
533 nrun = nrun + 1
534 110 CONTINUE
535 120 CONTINUE
536*
537* Print a summary of the results.
538*
539 CALL alasum( path, nout, nfail, nrun, nerrs )
540*
541 9999 FORMAT( ' N =', i5, ', type ', i2, ', test ', i2, ', ratio = ',
542 $ g12.5 )
543 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS =', i3,
544 $ ', type ', i2, ', test ', i2, ', ratio = ', g12.5 )
545 RETURN
546*
547* End of ZCHKPT
548*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
double precision function zlanht(norm, n, d, e)
ZLANHT returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanht.f:101
subroutine zlarnv(idist, iseed, n, x)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition zlarnv.f:99
subroutine dlarnv(idist, iseed, n, x)
DLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition dlarnv.f:97
subroutine zptcon(n, d, e, anorm, rcond, rwork, info)
ZPTCON
Definition zptcon.f:119
subroutine zptrfs(uplo, n, nrhs, d, e, df, ef, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPTRFS
Definition zptrfs.f:183
subroutine zpttrf(n, d, e, info)
ZPTTRF
Definition zpttrf.f:92
subroutine zpttrs(uplo, n, nrhs, d, e, b, ldb, info)
ZPTTRS
Definition zpttrs.f:121
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zerrgt(path, nunit)
ZERRGT
Definition zerrgt.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlaptm(uplo, n, nrhs, alpha, d, e, x, ldx, beta, b, ldb)
ZLAPTM
Definition zlaptm.f:129
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332
subroutine zptt01(n, d, e, df, ef, work, resid)
ZPTT01
Definition zptt01.f:92
subroutine zptt02(uplo, n, nrhs, d, e, x, ldx, b, ldb, resid)
ZPTT02
Definition zptt02.f:115
subroutine zptt05(n, nrhs, d, e, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPTT05
Definition zptt05.f:150
Here is the call graph for this function:
Here is the caller graph for this function: