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

ZCHKTR

Purpose:
 ZCHKTR tests ZTRTRI, -TRS, -RFS, and -CON, and ZLATRS
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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[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 maximum value of N in NVAL.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[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 165 of file zchktr.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: