LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zchkhe_rook ( 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( * )  AFAC,
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, dimension( * )  IWORK,
integer  NOUT 
)

ZCHKHE_ROOK

Purpose:
 ZCHKHE_ROOK tests ZHETRF_ROOK, -TRI_ROOK, -TRS_ROOK,
 and -CON_ROOK.
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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NBVAL)
          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 maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is CCOMPLEX*16 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is CCOMPLEX*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)
[out]IWORK
          IWORK is INTEGER array, dimension (2*NMAX)
[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
December 2016

Definition at line 174 of file zchkhe_rook.f.

174 *
175 * -- LAPACK test routine (version 3.7.0) --
176 * -- LAPACK is a software package provided by Univ. of Tennessee, --
177 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178 * December 2016
179 *
180 * .. Scalar Arguments ..
181  LOGICAL tsterr
182  INTEGER nmax, nn, nnb, nns, nout
183  DOUBLE PRECISION thresh
184 * ..
185 * .. Array Arguments ..
186  LOGICAL dotype( * )
187  INTEGER iwork( * ), nbval( * ), nsval( * ), nval( * )
188  DOUBLE PRECISION rwork( * )
189  COMPLEX*16 a( * ), afac( * ), ainv( * ), b( * ),
190  $ work( * ), x( * ), xact( * )
191 * ..
192 *
193 * =====================================================================
194 *
195 * .. Parameters ..
196  DOUBLE PRECISION zero, one
197  parameter ( zero = 0.0d+0, one = 1.0d+0 )
198  DOUBLE PRECISION onehalf
199  parameter ( onehalf = 0.5d+0 )
200  DOUBLE PRECISION eight, sevten
201  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
202  COMPLEX*16 czero
203  parameter ( czero = ( 0.0d+0, 0.0d+0 ) )
204  INTEGER ntypes
205  parameter ( ntypes = 10 )
206  INTEGER ntests
207  parameter ( ntests = 7 )
208 * ..
209 * .. Local Scalars ..
210  LOGICAL trfcon, zerot
211  CHARACTER dist, TYPE, uplo, xtype
212  CHARACTER*3 path, matpath
213  INTEGER i, i1, i2, imat, in, inb, info, ioff, irhs,
214  $ iuplo, izero, j, k, kl, ku, lda, lwork, mode,
215  $ n, nb, nerrs, nfail, nimat, nrhs, nrun, nt
216  DOUBLE PRECISION alpha, anorm, cndnum, const, sing_max,
217  $ sing_min, rcond, rcondc, dtemp
218 * ..
219 * .. Local Arrays ..
220  CHARACTER uplos( 2 )
221  INTEGER iseed( 4 ), iseedy( 4 )
222  DOUBLE PRECISION result( ntests )
223  COMPLEX*16 block( 2, 2 ), zdummy( 1 )
224 * ..
225 * .. External Functions ..
226  DOUBLE PRECISION zlange, zlanhe, dget06
227  EXTERNAL zlange, zlanhe, dget06
228 * ..
229 * .. External Subroutines ..
230  EXTERNAL alaerh, alahd, alasum, zerrhe, zgesvd, zget04,
234 * ..
235 * .. Intrinsic Functions ..
236  INTRINSIC conjg, max, min, sqrt
237 * ..
238 * .. Scalars in Common ..
239  LOGICAL lerr, ok
240  CHARACTER*32 srnamt
241  INTEGER infot, nunit
242 * ..
243 * .. Common blocks ..
244  COMMON / infoc / infot, nunit, ok, lerr
245  COMMON / srnamc / srnamt
246 * ..
247 * .. Data statements ..
248  DATA iseedy / 1988, 1989, 1990, 1991 /
249  DATA uplos / 'U', 'L' /
250 * ..
251 * .. Executable Statements ..
252 *
253 * Initialize constants and the random number seed.
254 *
255  alpha = ( one+sqrt( sevten ) ) / eight
256 *
257 * Test path
258 *
259  path( 1: 1 ) = 'Zomplex precision'
260  path( 2: 3 ) = 'HR'
261 *
262 * Path to generate matrices
263 *
264  matpath( 1: 1 ) = 'Zomplex precision'
265  matpath( 2: 3 ) = 'HE'
266 *
267  nrun = 0
268  nfail = 0
269  nerrs = 0
270  DO 10 i = 1, 4
271  iseed( i ) = iseedy( i )
272  10 CONTINUE
273 *
274 * Test the error exits
275 *
276  IF( tsterr )
277  $ CALL zerrhe( path, nout )
278  infot = 0
279 *
280 * Set the minimum block size for which the block routine should
281 * be used, which will be later returned by ILAENV
282 *
283  CALL xlaenv( 2, 2 )
284 *
285 * Do for each value of N in NVAL
286 *
287  DO 270 in = 1, nn
288  n = nval( in )
289  lda = max( n, 1 )
290  xtype = 'N'
291  nimat = ntypes
292  IF( n.LE.0 )
293  $ nimat = 1
294 *
295  izero = 0
296 *
297 * Do for each value of matrix type IMAT
298 *
299  DO 260 imat = 1, nimat
300 *
301 * Do the tests only if DOTYPE( IMAT ) is true.
302 *
303  IF( .NOT.dotype( imat ) )
304  $ GO TO 260
305 *
306 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
307 *
308  zerot = imat.GE.3 .AND. imat.LE.6
309  IF( zerot .AND. n.LT.imat-2 )
310  $ GO TO 260
311 *
312 * Do first for UPLO = 'U', then for UPLO = 'L'
313 *
314  DO 250 iuplo = 1, 2
315  uplo = uplos( iuplo )
316 *
317 * Begin generate the test matrix A.
318 *
319 * Set up parameters with ZLATB4 for the matrix generator
320 * based on the type of matrix to be generated.
321 *
322  CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
323  $ mode, cndnum, dist )
324 *
325 * Generate a matrix with ZLATMS.
326 *
327  srnamt = 'ZLATMS'
328  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
329  $ cndnum, anorm, kl, ku, uplo, a, lda,
330  $ work, info )
331 *
332 * Check error code from ZLATMS and handle error.
333 *
334  IF( info.NE.0 ) THEN
335  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
336  $ -1, -1, -1, imat, nfail, nerrs, nout )
337 *
338 * Skip all tests for this generated matrix
339 *
340  GO TO 250
341  END IF
342 *
343 * For matrix types 3-6, zero one or more rows and
344 * columns of the matrix to test that INFO is returned
345 * correctly.
346 *
347  IF( zerot ) THEN
348  IF( imat.EQ.3 ) THEN
349  izero = 1
350  ELSE IF( imat.EQ.4 ) THEN
351  izero = n
352  ELSE
353  izero = n / 2 + 1
354  END IF
355 *
356  IF( imat.LT.6 ) THEN
357 *
358 * Set row and column IZERO to zero.
359 *
360  IF( iuplo.EQ.1 ) THEN
361  ioff = ( izero-1 )*lda
362  DO 20 i = 1, izero - 1
363  a( ioff+i ) = czero
364  20 CONTINUE
365  ioff = ioff + izero
366  DO 30 i = izero, n
367  a( ioff ) = czero
368  ioff = ioff + lda
369  30 CONTINUE
370  ELSE
371  ioff = izero
372  DO 40 i = 1, izero - 1
373  a( ioff ) = czero
374  ioff = ioff + lda
375  40 CONTINUE
376  ioff = ioff - izero
377  DO 50 i = izero, n
378  a( ioff+i ) = czero
379  50 CONTINUE
380  END IF
381  ELSE
382  IF( iuplo.EQ.1 ) THEN
383 *
384 * Set the first IZERO rows and columns to zero.
385 *
386  ioff = 0
387  DO 70 j = 1, n
388  i2 = min( j, izero )
389  DO 60 i = 1, i2
390  a( ioff+i ) = czero
391  60 CONTINUE
392  ioff = ioff + lda
393  70 CONTINUE
394  ELSE
395 *
396 * Set the last IZERO rows and columns to zero.
397 *
398  ioff = 0
399  DO 90 j = 1, n
400  i1 = max( j, izero )
401  DO 80 i = i1, n
402  a( ioff+i ) = czero
403  80 CONTINUE
404  ioff = ioff + lda
405  90 CONTINUE
406  END IF
407  END IF
408  ELSE
409  izero = 0
410  END IF
411 *
412 * End generate the test matrix A.
413 *
414 *
415 * Do for each value of NB in NBVAL
416 *
417  DO 240 inb = 1, nnb
418 *
419 * Set the optimal blocksize, which will be later
420 * returned by ILAENV.
421 *
422  nb = nbval( inb )
423  CALL xlaenv( 1, nb )
424 *
425 * Copy the test matrix A into matrix AFAC which
426 * will be factorized in place. This is needed to
427 * preserve the test matrix A for subsequent tests.
428 *
429  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
430 *
431 * Compute the L*D*L**T or U*D*U**T factorization of the
432 * matrix. IWORK stores details of the interchanges and
433 * the block structure of D. AINV is a work array for
434 * block factorization, LWORK is the length of AINV.
435 *
436  lwork = max( 2, nb )*lda
437  srnamt = 'ZHETRF_ROOK'
438  CALL zhetrf_rook( uplo, n, afac, lda, iwork, ainv,
439  $ lwork, info )
440 *
441 * Adjust the expected value of INFO to account for
442 * pivoting.
443 *
444  k = izero
445  IF( k.GT.0 ) THEN
446  100 CONTINUE
447  IF( iwork( k ).LT.0 ) THEN
448  IF( iwork( k ).NE.-k ) THEN
449  k = -iwork( k )
450  GO TO 100
451  END IF
452  ELSE IF( iwork( k ).NE.k ) THEN
453  k = iwork( k )
454  GO TO 100
455  END IF
456  END IF
457 *
458 * Check error code from ZHETRF_ROOK and handle error.
459 *
460  IF( info.NE.k)
461  $ CALL alaerh( path, 'ZHETRF_ROOK', info, k,
462  $ uplo, n, n, -1, -1, nb, imat,
463  $ nfail, nerrs, nout )
464 *
465 * Set the condition estimate flag if the INFO is not 0.
466 *
467  IF( info.NE.0 ) THEN
468  trfcon = .true.
469  ELSE
470  trfcon = .false.
471  END IF
472 *
473 *+ TEST 1
474 * Reconstruct matrix from factors and compute residual.
475 *
476  CALL zhet01_rook( uplo, n, a, lda, afac, lda, iwork,
477  $ ainv, lda, rwork, result( 1 ) )
478  nt = 1
479 *
480 *+ TEST 2
481 * Form the inverse and compute the residual,
482 * if the factorization was competed without INFO > 0
483 * (i.e. there is no zero rows and columns).
484 * Do it only for the first block size.
485 *
486  IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
487  CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
488  srnamt = 'ZHETRI_ROOK'
489  CALL zhetri_rook( uplo, n, ainv, lda, iwork, work,
490  $ info )
491 *
492 * Check error code from ZHETRI_ROOK and handle error.
493 *
494  IF( info.NE.0 )
495  $ CALL alaerh( path, 'ZHETRI_ROOK', info, -1,
496  $ uplo, n, n, -1, -1, -1, imat,
497  $ nfail, nerrs, nout )
498 *
499 * Compute the residual for a Hermitian matrix times
500 * its inverse.
501 *
502  CALL zpot03( uplo, n, a, lda, ainv, lda, work, lda,
503  $ rwork, rcondc, result( 2 ) )
504  nt = 2
505  END IF
506 *
507 * Print information about the tests that did not pass
508 * the threshold.
509 *
510  DO 110 k = 1, nt
511  IF( result( k ).GE.thresh ) THEN
512  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
513  $ CALL alahd( nout, path )
514  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
515  $ result( k )
516  nfail = nfail + 1
517  END IF
518  110 CONTINUE
519  nrun = nrun + nt
520 *
521 *+ TEST 3
522 * Compute largest element in U or L
523 *
524  result( 3 ) = zero
525  dtemp = zero
526 *
527  const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) ) /
528  $ ( one-alpha )
529 *
530  IF( iuplo.EQ.1 ) THEN
531 *
532 * Compute largest element in U
533 *
534  k = n
535  120 CONTINUE
536  IF( k.LE.1 )
537  $ GO TO 130
538 *
539  IF( iwork( k ).GT.zero ) THEN
540 *
541 * Get max absolute value from elements
542 * in column k in U
543 *
544  dtemp = zlange( 'M', k-1, 1,
545  $ afac( ( k-1 )*lda+1 ), lda, rwork )
546  ELSE
547 *
548 * Get max absolute value from elements
549 * in columns k and k-1 in U
550 *
551  dtemp = zlange( 'M', k-2, 2,
552  $ afac( ( k-2 )*lda+1 ), lda, rwork )
553  k = k - 1
554 *
555  END IF
556 *
557 * DTEMP should be bounded by CONST
558 *
559  dtemp = dtemp - const + thresh
560  IF( dtemp.GT.result( 3 ) )
561  $ result( 3 ) = dtemp
562 *
563  k = k - 1
564 *
565  GO TO 120
566  130 CONTINUE
567 *
568  ELSE
569 *
570 * Compute largest element in L
571 *
572  k = 1
573  140 CONTINUE
574  IF( k.GE.n )
575  $ GO TO 150
576 *
577  IF( iwork( k ).GT.zero ) THEN
578 *
579 * Get max absolute value from elements
580 * in column k in L
581 *
582  dtemp = zlange( 'M', n-k, 1,
583  $ afac( ( k-1 )*lda+k+1 ), lda, rwork )
584  ELSE
585 *
586 * Get max absolute value from elements
587 * in columns k and k+1 in L
588 *
589  dtemp = zlange( 'M', n-k-1, 2,
590  $ afac( ( k-1 )*lda+k+2 ), lda, rwork )
591  k = k + 1
592 *
593  END IF
594 *
595 * DTEMP should be bounded by CONST
596 *
597  dtemp = dtemp - const + thresh
598  IF( dtemp.GT.result( 3 ) )
599  $ result( 3 ) = dtemp
600 *
601  k = k + 1
602 *
603  GO TO 140
604  150 CONTINUE
605  END IF
606 *
607 *
608 *+ TEST 4
609 * Compute largest 2-Norm (condition number)
610 * of 2-by-2 diag blocks
611 *
612  result( 4 ) = zero
613  dtemp = zero
614 *
615  const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) )*
616  $ ( ( one + alpha ) / ( one - alpha ) )
617  CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
618 *
619  IF( iuplo.EQ.1 ) THEN
620 *
621 * Loop backward for UPLO = 'U'
622 *
623  k = n
624  160 CONTINUE
625  IF( k.LE.1 )
626  $ GO TO 170
627 *
628  IF( iwork( k ).LT.zero ) THEN
629 *
630 * Get the two singular values
631 * (real and non-negative) of a 2-by-2 block,
632 * store them in RWORK array
633 *
634  block( 1, 1 ) = afac( ( k-2 )*lda+k-1 )
635  block( 1, 2 ) = afac( (k-1)*lda+k-1 )
636  block( 2, 1 ) = conjg( block( 1, 2 ) )
637  block( 2, 2 ) = afac( (k-1)*lda+k )
638 *
639  CALL zgesvd( 'N', 'N', 2, 2, block, 2, rwork,
640  $ zdummy, 1, zdummy, 1,
641  $ work, 6, rwork( 3 ), info )
642 *
643 *
644  sing_max = rwork( 1 )
645  sing_min = rwork( 2 )
646 *
647  dtemp = sing_max / sing_min
648 *
649 * DTEMP should be bounded by CONST
650 *
651  dtemp = dtemp - const + thresh
652  IF( dtemp.GT.result( 4 ) )
653  $ result( 4 ) = dtemp
654  k = k - 1
655 *
656  END IF
657 *
658  k = k - 1
659 *
660  GO TO 160
661  170 CONTINUE
662 *
663  ELSE
664 *
665 * Loop forward for UPLO = 'L'
666 *
667  k = 1
668  180 CONTINUE
669  IF( k.GE.n )
670  $ GO TO 190
671 *
672  IF( iwork( k ).LT.zero ) THEN
673 *
674 * Get the two singular values
675 * (real and non-negative) of a 2-by-2 block,
676 * store them in RWORK array
677 *
678  block( 1, 1 ) = afac( ( k-1 )*lda+k )
679  block( 2, 1 ) = afac( ( k-1 )*lda+k+1 )
680  block( 1, 2 ) = conjg( block( 2, 1 ) )
681  block( 2, 2 ) = afac( k*lda+k+1 )
682 *
683  CALL zgesvd( 'N', 'N', 2, 2, block, 2, rwork,
684  $ zdummy, 1, zdummy, 1,
685  $ work, 6, rwork(3), info )
686 *
687  sing_max = rwork( 1 )
688  sing_min = rwork( 2 )
689 *
690  dtemp = sing_max / sing_min
691 *
692 * DTEMP should be bounded by CONST
693 *
694  dtemp = dtemp - const + thresh
695  IF( dtemp.GT.result( 4 ) )
696  $ result( 4 ) = dtemp
697  k = k + 1
698 *
699  END IF
700 *
701  k = k + 1
702 *
703  GO TO 180
704  190 CONTINUE
705  END IF
706 *
707 * Print information about the tests that did not pass
708 * the threshold.
709 *
710  DO 200 k = 3, 4
711  IF( result( k ).GE.thresh ) THEN
712  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
713  $ CALL alahd( nout, path )
714  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
715  $ result( k )
716  nfail = nfail + 1
717  END IF
718  200 CONTINUE
719  nrun = nrun + 2
720 *
721 * Skip the other tests if this is not the first block
722 * size.
723 *
724  IF( inb.GT.1 )
725  $ GO TO 240
726 *
727 * Do only the condition estimate if INFO is not 0.
728 *
729  IF( trfcon ) THEN
730  rcondc = zero
731  GO TO 230
732  END IF
733 *
734 * Do for each value of NRHS in NSVAL.
735 *
736  DO 220 irhs = 1, nns
737  nrhs = nsval( irhs )
738 *
739 * Begin loop over NRHS values
740 *
741 *
742 *+ TEST 5 ( Using TRS_ROOK)
743 * Solve and compute residual for A * X = B.
744 *
745 * Choose a set of NRHS random solution vectors
746 * stored in XACT and set up the right hand side B
747 *
748  srnamt = 'ZLARHS'
749  CALL zlarhs( matpath, xtype, uplo, ' ', n, n,
750  $ kl, ku, nrhs, a, lda, xact, lda,
751  $ b, lda, iseed, info )
752  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
753 *
754  srnamt = 'ZHETRS_ROOK'
755  CALL zhetrs_rook( uplo, n, nrhs, afac, lda, iwork,
756  $ x, lda, info )
757 *
758 * Check error code from ZHETRS_ROOK and handle error.
759 *
760  IF( info.NE.0 )
761  $ CALL alaerh( path, 'ZHETRS_ROOK', info, 0,
762  $ uplo, n, n, -1, -1, nrhs, imat,
763  $ nfail, nerrs, nout )
764 *
765  CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
766 *
767 * Compute the residual for the solution
768 *
769  CALL zpot02( uplo, n, nrhs, a, lda, x, lda, work,
770  $ lda, rwork, result( 5 ) )
771 *
772 *+ TEST 6
773 * Check solution from generated exact solution.
774 *
775  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
776  $ result( 6 ) )
777 *
778 * Print information about the tests that did not pass
779 * the threshold.
780 *
781  DO 210 k = 5, 6
782  IF( result( k ).GE.thresh ) THEN
783  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
784  $ CALL alahd( nout, path )
785  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
786  $ imat, k, result( k )
787  nfail = nfail + 1
788  END IF
789  210 CONTINUE
790  nrun = nrun + 2
791 *
792 * End do for each value of NRHS in NSVAL.
793 *
794  220 CONTINUE
795 *
796 *+ TEST 7
797 * Get an estimate of RCOND = 1/CNDNUM.
798 *
799  230 CONTINUE
800  anorm = zlanhe( '1', uplo, n, a, lda, rwork )
801  srnamt = 'ZHECON_ROOK'
802  CALL zhecon_rook( uplo, n, afac, lda, iwork, anorm,
803  $ rcond, work, info )
804 *
805 * Check error code from ZHECON_ROOK and handle error.
806 *
807  IF( info.NE.0 )
808  $ CALL alaerh( path, 'ZHECON_ROOK', info, 0,
809  $ uplo, n, n, -1, -1, -1, imat,
810  $ nfail, nerrs, nout )
811 *
812 * Compute the test ratio to compare values of RCOND
813 *
814  result( 7 ) = dget06( rcond, rcondc )
815 *
816 * Print information about the tests that did not pass
817 * the threshold.
818 *
819  IF( result( 7 ).GE.thresh ) THEN
820  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
821  $ CALL alahd( nout, path )
822  WRITE( nout, fmt = 9997 )uplo, n, imat, 7,
823  $ result( 7 )
824  nfail = nfail + 1
825  END IF
826  nrun = nrun + 1
827  240 CONTINUE
828 *
829  250 CONTINUE
830  260 CONTINUE
831  270 CONTINUE
832 *
833 * Print a summary of the results.
834 *
835  CALL alasum( path, nout, nfail, nrun, nerrs )
836 *
837  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
838  $ i2, ', test ', i2, ', ratio =', g12.5 )
839  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
840  $ i2, ', test ', i2, ', ratio =', g12.5 )
841  9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
842  $ ', test ', i2, ', ratio =', g12.5 )
843  RETURN
844 *
845 * End of ZCHKHE_ROOK
846 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:107
subroutine zhet01_rook(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
ZHET01_ROOK
Definition: zhet01_rook.f:127
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
double precision function zlanhe(NORM, UPLO, N, A, LDA, WORK)
ZLANHE 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 matrix.
Definition: zlanhe.f:126
subroutine zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zerrhe(PATH, NUNIT)
ZERRHE
Definition: zerrhe.f:57
subroutine zpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZPOT02
Definition: zpot02.f:129
subroutine zhecon_rook(UPLO, N, A, LDA, IPIV, ANORM, RCOND, WORK, INFO)
ZHECON_ROOK estimates the reciprocal of the condition number fort HE matrices using factorization ob...
Definition: zhecon_rook.f:141
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 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 xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zpot03(UPLO, N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
ZPOT03
Definition: zpot03.f:128
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zhetrs_rook(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
ZHETRS_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using fac...
Definition: zhetrs_rook.f:138
subroutine zhetrf_rook(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
ZHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
Definition: zhetrf_rook.f:214
double precision function zlange(NORM, M, N, A, LDA, WORK)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlange.f:117
subroutine zgesvd(JOBU, JOBVT, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, INFO)
ZGESVD computes the singular value decomposition (SVD) for GE matrices
Definition: zgesvd.f:216
subroutine zhetri_rook(UPLO, N, A, LDA, IPIV, WORK, INFO)
ZHETRI_ROOK computes the inverse of HE matrix using the factorization obtained with the bounded Bunch...
Definition: zhetri_rook.f:130
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
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: