LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 149 of file zchkpt.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: