LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zchktp ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
complex*16, dimension( * )  AP,
complex*16, dimension( * )  AINVP,
complex*16, dimension( * )  B,
complex*16, dimension( * )  X,
complex*16, dimension( * )  XACT,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  NOUT 
)

ZCHKTP

Purpose:
 ZCHKTP tests ZTPTRI, -TRS, -RFS, and -CON, and ZLATPS
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 column 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.
[in]NMAX
          NMAX is INTEGER
          The leading dimension of the work arrays.  NMAX >= the
          maximumm value of N in NVAL.
[out]AP
          AP is COMPLEX*16 array, dimension (NMAX*(NMAX+1)/2)
[out]AINVP
          AINVP is COMPLEX*16 array, dimension (NMAX*(NMAX+1)/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 153 of file zchktp.f.

153 *
154 * -- LAPACK test routine (version 3.4.0) --
155 * -- LAPACK is a software package provided by Univ. of Tennessee, --
156 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157 * November 2011
158 *
159 * .. Scalar Arguments ..
160  LOGICAL tsterr
161  INTEGER nmax, nn, nns, nout
162  DOUBLE PRECISION thresh
163 * ..
164 * .. Array Arguments ..
165  LOGICAL dotype( * )
166  INTEGER nsval( * ), nval( * )
167  DOUBLE PRECISION rwork( * )
168  COMPLEX*16 ainvp( * ), ap( * ), b( * ), work( * ), x( * ),
169  $ xact( * )
170 * ..
171 *
172 * =====================================================================
173 *
174 * .. Parameters ..
175  INTEGER ntype1, ntypes
176  parameter ( ntype1 = 10, ntypes = 18 )
177  INTEGER ntests
178  parameter ( ntests = 9 )
179  INTEGER ntran
180  parameter ( ntran = 3 )
181  DOUBLE PRECISION one, zero
182  parameter ( one = 1.0d+0, zero = 0.0d+0 )
183 * ..
184 * .. Local Scalars ..
185  CHARACTER diag, norm, trans, uplo, xtype
186  CHARACTER*3 path
187  INTEGER i, idiag, imat, in, info, irhs, itran, iuplo,
188  $ k, lap, lda, n, nerrs, nfail, nrhs, nrun
189  DOUBLE PRECISION ainvnm, anorm, rcond, rcondc, rcondi, rcondo,
190  $ scale
191 * ..
192 * .. Local Arrays ..
193  CHARACTER transs( ntran ), uplos( 2 )
194  INTEGER iseed( 4 ), iseedy( 4 )
195  DOUBLE PRECISION result( ntests )
196 * ..
197 * .. External Functions ..
198  LOGICAL lsame
199  DOUBLE PRECISION zlantp
200  EXTERNAL lsame, zlantp
201 * ..
202 * .. External Subroutines ..
203  EXTERNAL alaerh, alahd, alasum, zcopy, zerrtr, zget04,
206  $ ztptrs
207 * ..
208 * .. Scalars in Common ..
209  LOGICAL lerr, ok
210  CHARACTER*32 srnamt
211  INTEGER infot, iounit
212 * ..
213 * .. Common blocks ..
214  COMMON / infoc / infot, iounit, ok, lerr
215  COMMON / srnamc / srnamt
216 * ..
217 * .. Intrinsic Functions ..
218  INTRINSIC max
219 * ..
220 * .. Data statements ..
221  DATA iseedy / 1988, 1989, 1990, 1991 /
222  DATA uplos / 'U', 'L' / , transs / 'N', 'T', 'C' /
223 * ..
224 * .. Executable Statements ..
225 *
226 * Initialize constants and the random number seed.
227 *
228  path( 1: 1 ) = 'Zomplex precision'
229  path( 2: 3 ) = 'TP'
230  nrun = 0
231  nfail = 0
232  nerrs = 0
233  DO 10 i = 1, 4
234  iseed( i ) = iseedy( i )
235  10 CONTINUE
236 *
237 * Test the error exits
238 *
239  IF( tsterr )
240  $ CALL zerrtr( path, nout )
241  infot = 0
242 *
243  DO 110 in = 1, nn
244 *
245 * Do for each value of N in NVAL
246 *
247  n = nval( in )
248  lda = max( 1, n )
249  lap = lda*( lda+1 ) / 2
250  xtype = 'N'
251 *
252  DO 70 imat = 1, ntype1
253 *
254 * Do the tests only if DOTYPE( IMAT ) is true.
255 *
256  IF( .NOT.dotype( imat ) )
257  $ GO TO 70
258 *
259  DO 60 iuplo = 1, 2
260 *
261 * Do first for UPLO = 'U', then for UPLO = 'L'
262 *
263  uplo = uplos( iuplo )
264 *
265 * Call ZLATTP to generate a triangular test matrix.
266 *
267  srnamt = 'ZLATTP'
268  CALL zlattp( imat, uplo, 'No transpose', diag, iseed, n,
269  $ ap, x, work, rwork, info )
270 *
271 * Set IDIAG = 1 for non-unit matrices, 2 for unit.
272 *
273  IF( lsame( diag, 'N' ) ) THEN
274  idiag = 1
275  ELSE
276  idiag = 2
277  END IF
278 *
279 *+ TEST 1
280 * Form the inverse of A.
281 *
282  IF( n.GT.0 )
283  $ CALL zcopy( lap, ap, 1, ainvp, 1 )
284  srnamt = 'ZTPTRI'
285  CALL ztptri( uplo, diag, n, ainvp, info )
286 *
287 * Check error code from ZTPTRI.
288 *
289  IF( info.NE.0 )
290  $ CALL alaerh( path, 'ZTPTRI', info, 0, uplo // diag, n,
291  $ n, -1, -1, -1, imat, nfail, nerrs, nout )
292 *
293 * Compute the infinity-norm condition number of A.
294 *
295  anorm = zlantp( 'I', uplo, diag, n, ap, rwork )
296  ainvnm = zlantp( 'I', uplo, diag, n, ainvp, rwork )
297  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
298  rcondi = one
299  ELSE
300  rcondi = ( one / anorm ) / ainvnm
301  END IF
302 *
303 * Compute the residual for the triangular matrix times its
304 * inverse. Also compute the 1-norm condition number of A.
305 *
306  CALL ztpt01( uplo, diag, n, ap, ainvp, rcondo, rwork,
307  $ result( 1 ) )
308 *
309 * Print the test ratio if it is .GE. THRESH.
310 *
311  IF( result( 1 ).GE.thresh ) THEN
312  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
313  $ CALL alahd( nout, path )
314  WRITE( nout, fmt = 9999 )uplo, diag, n, imat, 1,
315  $ result( 1 )
316  nfail = nfail + 1
317  END IF
318  nrun = nrun + 1
319 *
320  DO 40 irhs = 1, nns
321  nrhs = nsval( irhs )
322  xtype = 'N'
323 *
324  DO 30 itran = 1, ntran
325 *
326 * Do for op(A) = A, A**T, or A**H.
327 *
328  trans = transs( itran )
329  IF( itran.EQ.1 ) THEN
330  norm = 'O'
331  rcondc = rcondo
332  ELSE
333  norm = 'I'
334  rcondc = rcondi
335  END IF
336 *
337 *+ TEST 2
338 * Solve and compute residual for op(A)*x = b.
339 *
340  srnamt = 'ZLARHS'
341  CALL zlarhs( path, xtype, uplo, trans, n, n, 0,
342  $ idiag, nrhs, ap, lap, xact, lda, b,
343  $ lda, iseed, info )
344  xtype = 'C'
345  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
346 *
347  srnamt = 'ZTPTRS'
348  CALL ztptrs( uplo, trans, diag, n, nrhs, ap, x,
349  $ lda, info )
350 *
351 * Check error code from ZTPTRS.
352 *
353  IF( info.NE.0 )
354  $ CALL alaerh( path, 'ZTPTRS', info, 0,
355  $ uplo // trans // diag, n, n, -1,
356  $ -1, -1, imat, nfail, nerrs, nout )
357 *
358  CALL ztpt02( uplo, trans, diag, n, nrhs, ap, x,
359  $ lda, b, lda, work, rwork,
360  $ result( 2 ) )
361 *
362 *+ TEST 3
363 * Check solution from generated exact solution.
364 *
365  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
366  $ result( 3 ) )
367 *
368 *+ TESTS 4, 5, and 6
369 * Use iterative refinement to improve the solution and
370 * compute error bounds.
371 *
372  srnamt = 'ZTPRFS'
373  CALL ztprfs( uplo, trans, diag, n, nrhs, ap, b,
374  $ lda, x, lda, rwork, rwork( nrhs+1 ),
375  $ work, rwork( 2*nrhs+1 ), info )
376 *
377 * Check error code from ZTPRFS.
378 *
379  IF( info.NE.0 )
380  $ CALL alaerh( path, 'ZTPRFS', info, 0,
381  $ uplo // trans // diag, n, n, -1,
382  $ -1, nrhs, imat, nfail, nerrs,
383  $ nout )
384 *
385  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
386  $ result( 4 ) )
387  CALL ztpt05( uplo, trans, diag, n, nrhs, ap, b,
388  $ lda, x, lda, xact, lda, rwork,
389  $ rwork( nrhs+1 ), result( 5 ) )
390 *
391 * Print information about the tests that did not pass
392 * the threshold.
393 *
394  DO 20 k = 2, 6
395  IF( result( k ).GE.thresh ) THEN
396  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
397  $ CALL alahd( nout, path )
398  WRITE( nout, fmt = 9998 )uplo, trans, diag,
399  $ n, nrhs, imat, k, result( k )
400  nfail = nfail + 1
401  END IF
402  20 CONTINUE
403  nrun = nrun + 5
404  30 CONTINUE
405  40 CONTINUE
406 *
407 *+ TEST 7
408 * Get an estimate of RCOND = 1/CNDNUM.
409 *
410  DO 50 itran = 1, 2
411  IF( itran.EQ.1 ) THEN
412  norm = 'O'
413  rcondc = rcondo
414  ELSE
415  norm = 'I'
416  rcondc = rcondi
417  END IF
418  srnamt = 'ZTPCON'
419  CALL ztpcon( norm, uplo, diag, n, ap, rcond, work,
420  $ rwork, info )
421 *
422 * Check error code from ZTPCON.
423 *
424  IF( info.NE.0 )
425  $ CALL alaerh( path, 'ZTPCON', info, 0,
426  $ norm // uplo // diag, n, n, -1, -1,
427  $ -1, imat, nfail, nerrs, nout )
428 *
429  CALL ztpt06( rcond, rcondc, uplo, diag, n, ap, rwork,
430  $ result( 7 ) )
431 *
432 * Print the test ratio if it is .GE. THRESH.
433 *
434  IF( result( 7 ).GE.thresh ) THEN
435  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
436  $ CALL alahd( nout, path )
437  WRITE( nout, fmt = 9997 ) 'ZTPCON', norm, uplo,
438  $ diag, n, imat, 7, result( 7 )
439  nfail = nfail + 1
440  END IF
441  nrun = nrun + 1
442  50 CONTINUE
443  60 CONTINUE
444  70 CONTINUE
445 *
446 * Use pathological test matrices to test ZLATPS.
447 *
448  DO 100 imat = ntype1 + 1, ntypes
449 *
450 * Do the tests only if DOTYPE( IMAT ) is true.
451 *
452  IF( .NOT.dotype( imat ) )
453  $ GO TO 100
454 *
455  DO 90 iuplo = 1, 2
456 *
457 * Do first for UPLO = 'U', then for UPLO = 'L'
458 *
459  uplo = uplos( iuplo )
460  DO 80 itran = 1, ntran
461 *
462 * Do for op(A) = A, A**T, or A**H.
463 *
464  trans = transs( itran )
465 *
466 * Call ZLATTP to generate a triangular test matrix.
467 *
468  srnamt = 'ZLATTP'
469  CALL zlattp( imat, uplo, trans, diag, iseed, n, ap, x,
470  $ work, rwork, info )
471 *
472 *+ TEST 8
473 * Solve the system op(A)*x = b.
474 *
475  srnamt = 'ZLATPS'
476  CALL zcopy( n, x, 1, b, 1 )
477  CALL zlatps( uplo, trans, diag, 'N', n, ap, b, scale,
478  $ rwork, info )
479 *
480 * Check error code from ZLATPS.
481 *
482  IF( info.NE.0 )
483  $ CALL alaerh( path, 'ZLATPS', info, 0,
484  $ uplo // trans // diag // 'N', n, n,
485  $ -1, -1, -1, imat, nfail, nerrs, nout )
486 *
487  CALL ztpt03( uplo, trans, diag, n, 1, ap, scale,
488  $ rwork, one, b, lda, x, lda, work,
489  $ result( 8 ) )
490 *
491 *+ TEST 9
492 * Solve op(A)*x = b again with NORMIN = 'Y'.
493 *
494  CALL zcopy( n, x, 1, b( n+1 ), 1 )
495  CALL zlatps( uplo, trans, diag, 'Y', n, ap, b( n+1 ),
496  $ scale, rwork, info )
497 *
498 * Check error code from ZLATPS.
499 *
500  IF( info.NE.0 )
501  $ CALL alaerh( path, 'ZLATPS', info, 0,
502  $ uplo // trans // diag // 'Y', n, n,
503  $ -1, -1, -1, imat, nfail, nerrs, nout )
504 *
505  CALL ztpt03( uplo, trans, diag, n, 1, ap, scale,
506  $ rwork, one, b( n+1 ), lda, x, lda, work,
507  $ result( 9 ) )
508 *
509 * Print information about the tests that did not pass
510 * the threshold.
511 *
512  IF( result( 8 ).GE.thresh ) THEN
513  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
514  $ CALL alahd( nout, path )
515  WRITE( nout, fmt = 9996 )'ZLATPS', uplo, trans,
516  $ diag, 'N', n, imat, 8, result( 8 )
517  nfail = nfail + 1
518  END IF
519  IF( result( 9 ).GE.thresh ) THEN
520  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
521  $ CALL alahd( nout, path )
522  WRITE( nout, fmt = 9996 )'ZLATPS', uplo, trans,
523  $ diag, 'Y', n, imat, 9, result( 9 )
524  nfail = nfail + 1
525  END IF
526  nrun = nrun + 2
527  80 CONTINUE
528  90 CONTINUE
529  100 CONTINUE
530  110 CONTINUE
531 *
532 * Print a summary of the results.
533 *
534  CALL alasum( path, nout, nfail, nrun, nerrs )
535 *
536  9999 FORMAT( ' UPLO=''', a1, ''', DIAG=''', a1, ''', N=', i5,
537  $ ', type ', i2, ', test(', i2, ')= ', g12.5 )
538  9998 FORMAT( ' UPLO=''', a1, ''', TRANS=''', a1, ''', DIAG=''', a1,
539  $ ''', N=', i5, ''', NRHS=', i5, ', type ', i2, ', test(',
540  $ i2, ')= ', g12.5 )
541  9997 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''',',
542  $ i5, ', ... ), type ', i2, ', test(', i2, ')=', g12.5 )
543  9996 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''', ''',
544  $ a1, ''',', i5, ', ... ), type ', i2, ', test(', i2, ')=',
545  $ g12.5 )
546  RETURN
547 *
548 * End of ZCHKTP
549 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine zerrtr(PATH, NUNIT)
ZERRTR
Definition: zerrtr.f:56
subroutine ztptri(UPLO, DIAG, N, AP, INFO)
ZTPTRI
Definition: ztptri.f:119
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 zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine ztpt02(UPLO, TRANS, DIAG, N, NRHS, AP, X, LDX, B, LDB, WORK, RWORK, RESID)
ZTPT02
Definition: ztpt02.f:151
subroutine ztpt05(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
ZTPT05
Definition: ztpt05.f:177
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
subroutine zlattp(IMAT, UPLO, TRANS, DIAG, ISEED, N, AP, B, WORK, RWORK, INFO)
ZLATTP
Definition: zlattp.f:133
double precision function zlantp(NORM, UPLO, DIAG, N, AP, WORK)
ZLANTP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a triangular matrix supplied in packed form.
Definition: zlantp.f:127
subroutine ztpcon(NORM, UPLO, DIAG, N, AP, RCOND, WORK, RWORK, INFO)
ZTPCON
Definition: ztpcon.f:132
subroutine ztpt03(UPLO, TRANS, DIAG, N, NRHS, AP, SCALE, CNORM, TSCAL, X, LDX, B, LDB, WORK, RESID)
ZTPT03
Definition: ztpt03.f:164
subroutine ztptrs(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, INFO)
ZTPTRS
Definition: ztptrs.f:132
subroutine ztpt01(UPLO, DIAG, N, AP, AINVP, RCOND, RWORK, RESID)
ZTPT01
Definition: ztpt01.f:111
subroutine zlatps(UPLO, TRANS, DIAG, NORMIN, N, AP, X, SCALE, CNORM, INFO)
ZLATPS solves a triangular system of equations with the matrix held in packed storage.
Definition: zlatps.f:233
subroutine ztprfs(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
ZTPRFS
Definition: ztprfs.f:176
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75
subroutine ztpt06(RCOND, RCONDC, UPLO, DIAG, N, AP, RWORK, RAT)
ZTPT06
Definition: ztpt06.f:114

Here is the call graph for this function:

Here is the caller graph for this function: