LAPACK  3.6.1
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
November 2015

Definition at line 174 of file zchkhe_rook.f.

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