LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine cdrvst2stg ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
real  THRESH,
integer  NOUNIT,
complex, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  D1,
real, dimension( * )  D2,
real, dimension( * )  D3,
real, dimension( * )  WA1,
real, dimension( * )  WA2,
real, dimension( * )  WA3,
complex, dimension( ldu, * )  U,
integer  LDU,
complex, dimension( ldu, * )  V,
complex, dimension( * )  TAU,
complex, dimension( ldu, * )  Z,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
real, dimension( * )  RESULT,
integer  INFO 
)

CDRVST2STG

Purpose:
      CDRVST2STG  checks the Hermitian eigenvalue problem drivers.

              CHEEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix,
              using a divide-and-conquer algorithm.

              CHEEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix.

              CHEEVR computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix
              using the Relatively Robust Representation where it can.

              CHPEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage, using a divide-and-conquer algorithm.

              CHPEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage.

              CHBEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix,
              using a divide-and-conquer algorithm.

              CHBEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix.

              CHEEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix.

              CHPEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage.

              CHBEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix.

      When CDRVST2STG is called, a number of matrix "sizes" ("n's") and a
      number of matrix "types" are specified.  For each size ("n")
      and each type of matrix, one matrix will be generated and used
      to test the appropriate drivers.  For each matrix and each
      driver routine called, the following tests will be performed:

      (1)     | A - Z D Z' | / ( |A| n ulp )

      (2)     | I - Z Z' | / ( n ulp )

      (3)     | D1 - D2 | / ( |D1| ulp )

      where Z is the matrix of eigenvectors returned when the
      eigenvector option is given and D1 and D2 are the eigenvalues
      returned with and without the eigenvector option.

      The "sizes" are specified by an array NN(1:NSIZES); the value of
      each element NN(j) specifies one size.
      The "types" are specified by a logical array DOTYPE( 1:NTYPES );
      if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
      Currently, the list of possible types is:

      (1)  The zero matrix.
      (2)  The identity matrix.

      (3)  A diagonal matrix with evenly spaced entries
           1, ..., ULP  and random signs.
           (ULP = (first number larger than 1) - 1 )
      (4)  A diagonal matrix with geometrically spaced entries
           1, ..., ULP  and random signs.
      (5)  A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
           and random signs.

      (6)  Same as (4), but multiplied by SQRT( overflow threshold )
      (7)  Same as (4), but multiplied by SQRT( underflow threshold )

      (8)  A matrix of the form  U* D U, where U is unitary and
           D has evenly spaced entries 1, ..., ULP with random signs
           on the diagonal.

      (9)  A matrix of the form  U* D U, where U is unitary and
           D has geometrically spaced entries 1, ..., ULP with random
           signs on the diagonal.

      (10) A matrix of the form  U* D U, where U is unitary and
           D has "clustered" entries 1, ULP,..., ULP with random
           signs on the diagonal.

      (11) Same as (8), but multiplied by SQRT( overflow threshold )
      (12) Same as (8), but multiplied by SQRT( underflow threshold )

      (13) Symmetric matrix with random entries chosen from (-1,1).
      (14) Same as (13), but multiplied by SQRT( overflow threshold )
      (15) Same as (13), but multiplied by SQRT( underflow threshold )
      (16) A band matrix with half bandwidth randomly chosen between
           0 and N-1, with evenly spaced eigenvalues 1, ..., ULP
           with random signs.
      (17) Same as (16), but multiplied by SQRT( overflow threshold )
      (18) Same as (16), but multiplied by SQRT( underflow threshold )
  NSIZES  INTEGER
          The number of sizes of matrices to use.  If it is zero,
          CDRVST2STG does nothing.  It must be at least zero.
          Not modified.

  NN      INTEGER array, dimension (NSIZES)
          An array containing the sizes to be used for the matrices.
          Zero values will be skipped.  The values must be at least
          zero.
          Not modified.

  NTYPES  INTEGER
          The number of elements in DOTYPE.   If it is zero, CDRVST2STG
          does nothing.  It must be at least zero.  If it is MAXTYP+1
          and NSIZES is 1, then an additional type, MAXTYP+1 is
          defined, which is to use whatever matrix is in A.  This
          is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
          DOTYPE(MAXTYP+1) is .TRUE. .
          Not modified.

  DOTYPE  LOGICAL array, dimension (NTYPES)
          If DOTYPE(j) is .TRUE., then for each size in NN a
          matrix of that size and of type j will be generated.
          If NTYPES is smaller than the maximum number of types
          defined (PARAMETER MAXTYP), then types NTYPES+1 through
          MAXTYP will not be generated.  If NTYPES is larger
          than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
          will be ignored.
          Not modified.

  ISEED   INTEGER array, dimension (4)
          On entry ISEED specifies the seed of the random number
          generator. The array elements should be between 0 and 4095;
          if not they will be reduced mod 4096.  Also, ISEED(4) must
          be odd.  The random number generator uses a linear
          congruential sequence limited to small integers, and so
          should produce machine independent random numbers. The
          values of ISEED are changed on exit, and can be used in the
          next call to CDRVST2STG to continue the same random number
          sequence.
          Modified.

  THRESH  REAL
          A test will count as "failed" if the "error", computed as
          described above, exceeds THRESH.  Note that the error
          is scaled to be O(1), so THRESH should be a reasonably
          small multiple of 1, e.g., 10 or 100.  In particular,
          it should not depend on the precision (single vs. double)
          or the size of the matrix.  It must be at least zero.
          Not modified.

  NOUNIT  INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
          Not modified.

  A       COMPLEX array, dimension (LDA , max(NN))
          Used to hold the matrix whose eigenvalues are to be
          computed.  On exit, A contains the last matrix actually
          used.
          Modified.

  LDA     INTEGER
          The leading dimension of A.  It must be at
          least 1 and at least max( NN ).
          Not modified.

  D1      REAL             array, dimension (max(NN))
          The eigenvalues of A, as computed by CSTEQR simlutaneously
          with Z.  On exit, the eigenvalues in D1 correspond with the
          matrix in A.
          Modified.

  D2      REAL             array, dimension (max(NN))
          The eigenvalues of A, as computed by CSTEQR if Z is not
          computed.  On exit, the eigenvalues in D2 correspond with
          the matrix in A.
          Modified.

  D3      REAL             array, dimension (max(NN))
          The eigenvalues of A, as computed by SSTERF.  On exit, the
          eigenvalues in D3 correspond with the matrix in A.
          Modified.

  WA1     REAL array, dimension

  WA2     REAL array, dimension

  WA3     REAL array, dimension

  U       COMPLEX array, dimension (LDU, max(NN))
          The unitary matrix computed by CHETRD + CUNGC3.
          Modified.

  LDU     INTEGER
          The leading dimension of U, Z, and V.  It must be at
          least 1 and at least max( NN ).
          Not modified.

  V       COMPLEX array, dimension (LDU, max(NN))
          The Housholder vectors computed by CHETRD in reducing A to
          tridiagonal form.
          Modified.

  TAU     COMPLEX array, dimension (max(NN))
          The Householder factors computed by CHETRD in reducing A
          to tridiagonal form.
          Modified.

  Z       COMPLEX array, dimension (LDU, max(NN))
          The unitary matrix of eigenvectors computed by CHEEVD,
          CHEEVX, CHPEVD, CHPEVX, CHBEVD, and CHBEVX.
          Modified.

  WORK  - COMPLEX array of dimension ( LWORK )
           Workspace.
           Modified.

  LWORK - INTEGER
           The number of entries in WORK.  This must be at least
           2*max( NN(j), 2 )**2.
           Not modified.

  RWORK   REAL array, dimension (3*max(NN))
           Workspace.
           Modified.

  LRWORK - INTEGER
           The number of entries in RWORK.

  IWORK   INTEGER array, dimension (6*max(NN))
          Workspace.
          Modified.

  LIWORK - INTEGER
           The number of entries in IWORK.

  RESULT  REAL array, dimension (??)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid
          overflow.
          Modified.

  INFO    INTEGER
          If 0, then everything ran OK.
           -1: NSIZES < 0
           -2: Some NN(j) < 0
           -3: NTYPES < 0
           -5: THRESH < 0
           -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
          -16: LDU < 1 or LDU < NMAX.
          -21: LWORK too small.
          If  SLATMR, SLATMS, CHETRD, SORGC3, CSTEQR, SSTERF,
              or SORMC2 returns an error code, the
              absolute value of it is returned.
          Modified.

-----------------------------------------------------------------------

       Some Local Variables and Parameters:
       ---- ----- --------- --- ----------
       ZERO, ONE       Real 0 and 1.
       MAXTYP          The number of types defined.
       NTEST           The number of tests performed, or which can
                       be performed so far, for the current matrix.
       NTESTT          The total number of tests performed so far.
       NMAX            Largest value in NN.
       NMATS           The number of matrices generated so far.
       NERRS           The number of tests which have exceeded THRESH
                       so far (computed by SLAFTS).
       COND, IMODE     Values to be passed to the matrix generators.
       ANORM           Norm of A; passed to matrix generators.

       OVFL, UNFL      Overflow and underflow thresholds.
       ULP, ULPINV     Finest relative precision and its inverse.
       RTOVFL, RTUNFL  Square roots of the previous 2 values.
               The following four arrays decode JTYPE:
       KTYPE(j)        The general type (1-10) for type "j".
       KMODE(j)        The MODE value to be passed to the matrix
                       generator for type "j".
       KMAGN(j)        The order of magnitude ( O(1),
                       O(overflow^(1/2) ), O(underflow^(1/2) )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016

Definition at line 340 of file cdrvst2stg.f.

340 *
341 * -- LAPACK test routine (version 3.7.0) --
342 * -- LAPACK is a software package provided by Univ. of Tennessee, --
343 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
344 * December 2016
345 *
346 * .. Scalar Arguments ..
347  INTEGER info, lda, ldu, liwork, lrwork, lwork, nounit,
348  $ nsizes, ntypes
349  REAL thresh
350 * ..
351 * .. Array Arguments ..
352  LOGICAL dotype( * )
353  INTEGER iseed( 4 ), iwork( * ), nn( * )
354  REAL d1( * ), d2( * ), d3( * ), result( * ),
355  $ rwork( * ), wa1( * ), wa2( * ), wa3( * )
356  COMPLEX a( lda, * ), tau( * ), u( ldu, * ),
357  $ v( ldu, * ), work( * ), z( ldu, * )
358 * ..
359 *
360 * =====================================================================
361 *
362 *
363 * .. Parameters ..
364  REAL zero, one, two, ten
365  parameter ( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0,
366  $ ten = 10.0e+0 )
367  REAL half
368  parameter ( half = one / two )
369  COMPLEX czero, cone
370  parameter ( czero = ( 0.0e+0, 0.0e+0 ),
371  $ cone = ( 1.0e+0, 0.0e+0 ) )
372  INTEGER maxtyp
373  parameter ( maxtyp = 18 )
374 * ..
375 * .. Local Scalars ..
376  LOGICAL badnn
377  CHARACTER uplo
378  INTEGER i, idiag, ihbw, iinfo, il, imode, indwrk, indx,
379  $ irow, itemp, itype, iu, iuplo, j, j1, j2, jcol,
380  $ jsize, jtype, kd, lgn, liwedc, lrwedc, lwedc,
381  $ m, m2, m3, mtypes, n, nerrs, nmats, nmax,
382  $ ntest, ntestt
383  REAL abstol, aninv, anorm, cond, ovfl, rtovfl,
384  $ rtunfl, temp1, temp2, temp3, ulp, ulpinv, unfl,
385  $ vl, vu
386 * ..
387 * .. Local Arrays ..
388  INTEGER idumma( 1 ), ioldsd( 4 ), iseed2( 4 ),
389  $ iseed3( 4 ), kmagn( maxtyp ), kmode( maxtyp ),
390  $ ktype( maxtyp )
391 * ..
392 * .. External Functions ..
393  REAL slamch, slarnd, ssxt1
394  EXTERNAL slamch, slarnd, ssxt1
395 * ..
396 * .. External Subroutines ..
397  EXTERNAL alasvm, slabad, slafts, xerbla, chbev, chbevd,
402  $ chbevx_2stage, chetrd_2stage, chetrd_sy2sb,
403  $ chetrd_sb2st, clatmr, clatms
404 * ..
405 * .. Intrinsic Functions ..
406  INTRINSIC abs, REAL, int, log, max, min, sqrt
407 * ..
408 * .. Data statements ..
409  DATA ktype / 1, 2, 5*4, 5*5, 3*8, 3*9 /
410  DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
411  $ 2, 3, 1, 2, 3 /
412  DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
413  $ 0, 0, 4, 4, 4 /
414 * ..
415 * .. Executable Statements ..
416 *
417 * 1) Check for errors
418 *
419  ntestt = 0
420  info = 0
421 *
422  badnn = .false.
423  nmax = 1
424  DO 10 j = 1, nsizes
425  nmax = max( nmax, nn( j ) )
426  IF( nn( j ).LT.0 )
427  $ badnn = .true.
428  10 CONTINUE
429 *
430 * Check for errors
431 *
432  IF( nsizes.LT.0 ) THEN
433  info = -1
434  ELSE IF( badnn ) THEN
435  info = -2
436  ELSE IF( ntypes.LT.0 ) THEN
437  info = -3
438  ELSE IF( lda.LT.nmax ) THEN
439  info = -9
440  ELSE IF( ldu.LT.nmax ) THEN
441  info = -16
442  ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
443  info = -22
444  END IF
445 *
446  IF( info.NE.0 ) THEN
447  CALL xerbla( 'CDRVST2STG', -info )
448  RETURN
449  END IF
450 *
451 * Quick return if nothing to do
452 *
453  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
454  $ RETURN
455 *
456 * More Important constants
457 *
458  unfl = slamch( 'Safe minimum' )
459  ovfl = slamch( 'Overflow' )
460  CALL slabad( unfl, ovfl )
461  ulp = slamch( 'Epsilon' )*slamch( 'Base' )
462  ulpinv = one / ulp
463  rtunfl = sqrt( unfl )
464  rtovfl = sqrt( ovfl )
465 *
466 * Loop over sizes, types
467 *
468  DO 20 i = 1, 4
469  iseed2( i ) = iseed( i )
470  iseed3( i ) = iseed( i )
471  20 CONTINUE
472 *
473  nerrs = 0
474  nmats = 0
475 *
476  DO 1220 jsize = 1, nsizes
477  n = nn( jsize )
478  IF( n.GT.0 ) THEN
479  lgn = int( log( REAL( N ) ) / log( two ) )
480  IF( 2**lgn.LT.n )
481  $ lgn = lgn + 1
482  IF( 2**lgn.LT.n )
483  $ lgn = lgn + 1
484  lwedc = max( 2*n+n*n, 2*n*n )
485  lrwedc = 1 + 4*n + 2*n*lgn + 3*n**2
486  liwedc = 3 + 5*n
487  ELSE
488  lwedc = 2
489  lrwedc = 8
490  liwedc = 8
491  END IF
492  aninv = one / REAL( MAX( 1, N ) )
493 *
494  IF( nsizes.NE.1 ) THEN
495  mtypes = min( maxtyp, ntypes )
496  ELSE
497  mtypes = min( maxtyp+1, ntypes )
498  END IF
499 *
500  DO 1210 jtype = 1, mtypes
501  IF( .NOT.dotype( jtype ) )
502  $ GO TO 1210
503  nmats = nmats + 1
504  ntest = 0
505 *
506  DO 30 j = 1, 4
507  ioldsd( j ) = iseed( j )
508  30 CONTINUE
509 *
510 * 2) Compute "A"
511 *
512 * Control parameters:
513 *
514 * KMAGN KMODE KTYPE
515 * =1 O(1) clustered 1 zero
516 * =2 large clustered 2 identity
517 * =3 small exponential (none)
518 * =4 arithmetic diagonal, (w/ eigenvalues)
519 * =5 random log Hermitian, w/ eigenvalues
520 * =6 random (none)
521 * =7 random diagonal
522 * =8 random Hermitian
523 * =9 band Hermitian, w/ eigenvalues
524 *
525  IF( mtypes.GT.maxtyp )
526  $ GO TO 110
527 *
528  itype = ktype( jtype )
529  imode = kmode( jtype )
530 *
531 * Compute norm
532 *
533  GO TO ( 40, 50, 60 )kmagn( jtype )
534 *
535  40 CONTINUE
536  anorm = one
537  GO TO 70
538 *
539  50 CONTINUE
540  anorm = ( rtovfl*ulp )*aninv
541  GO TO 70
542 *
543  60 CONTINUE
544  anorm = rtunfl*n*ulpinv
545  GO TO 70
546 *
547  70 CONTINUE
548 *
549  CALL claset( 'Full', lda, n, czero, czero, a, lda )
550  iinfo = 0
551  cond = ulpinv
552 *
553 * Special Matrices -- Identity & Jordan block
554 *
555 * Zero
556 *
557  IF( itype.EQ.1 ) THEN
558  iinfo = 0
559 *
560  ELSE IF( itype.EQ.2 ) THEN
561 *
562 * Identity
563 *
564  DO 80 jcol = 1, n
565  a( jcol, jcol ) = anorm
566  80 CONTINUE
567 *
568  ELSE IF( itype.EQ.4 ) THEN
569 *
570 * Diagonal Matrix, [Eigen]values Specified
571 *
572  CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
573  $ anorm, 0, 0, 'N', a, lda, work, iinfo )
574 *
575  ELSE IF( itype.EQ.5 ) THEN
576 *
577 * Hermitian, eigenvalues specified
578 *
579  CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
580  $ anorm, n, n, 'N', a, lda, work, iinfo )
581 *
582  ELSE IF( itype.EQ.7 ) THEN
583 *
584 * Diagonal, random eigenvalues
585 *
586  CALL clatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
587  $ 'T', 'N', work( n+1 ), 1, one,
588  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
589  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
590 *
591  ELSE IF( itype.EQ.8 ) THEN
592 *
593 * Hermitian, random eigenvalues
594 *
595  CALL clatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
596  $ 'T', 'N', work( n+1 ), 1, one,
597  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
598  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
599 *
600  ELSE IF( itype.EQ.9 ) THEN
601 *
602 * Hermitian banded, eigenvalues specified
603 *
604  ihbw = int( ( n-1 )*slarnd( 1, iseed3 ) )
605  CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
606  $ anorm, ihbw, ihbw, 'Z', u, ldu, work,
607  $ iinfo )
608 *
609 * Store as dense matrix for most routines.
610 *
611  CALL claset( 'Full', lda, n, czero, czero, a, lda )
612  DO 100 idiag = -ihbw, ihbw
613  irow = ihbw - idiag + 1
614  j1 = max( 1, idiag+1 )
615  j2 = min( n, n+idiag )
616  DO 90 j = j1, j2
617  i = j - idiag
618  a( i, j ) = u( irow, j )
619  90 CONTINUE
620  100 CONTINUE
621  ELSE
622  iinfo = 1
623  END IF
624 *
625  IF( iinfo.NE.0 ) THEN
626  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
627  $ ioldsd
628  info = abs( iinfo )
629  RETURN
630  END IF
631 *
632  110 CONTINUE
633 *
634  abstol = unfl + unfl
635  IF( n.LE.1 ) THEN
636  il = 1
637  iu = n
638  ELSE
639  il = 1 + int( ( n-1 )*slarnd( 1, iseed2 ) )
640  iu = 1 + int( ( n-1 )*slarnd( 1, iseed2 ) )
641  IF( il.GT.iu ) THEN
642  itemp = il
643  il = iu
644  iu = itemp
645  END IF
646  END IF
647 *
648 * Perform tests storing upper or lower triangular
649 * part of matrix.
650 *
651  DO 1200 iuplo = 0, 1
652  IF( iuplo.EQ.0 ) THEN
653  uplo = 'L'
654  ELSE
655  uplo = 'U'
656  END IF
657 *
658 * Call CHEEVD and CHEEVX.
659 *
660  CALL clacpy( ' ', n, n, a, lda, v, ldu )
661 *
662  ntest = ntest + 1
663  CALL cheevd( 'V', uplo, n, a, ldu, d1, work, lwedc,
664  $ rwork, lrwedc, iwork, liwedc, iinfo )
665  IF( iinfo.NE.0 ) THEN
666  WRITE( nounit, fmt = 9999 )'CHEEVD(V,' // uplo //
667  $ ')', iinfo, n, jtype, ioldsd
668  info = abs( iinfo )
669  IF( iinfo.LT.0 ) THEN
670  RETURN
671  ELSE
672  result( ntest ) = ulpinv
673  result( ntest+1 ) = ulpinv
674  result( ntest+2 ) = ulpinv
675  GO TO 130
676  END IF
677  END IF
678 *
679 * Do tests 1 and 2.
680 *
681  CALL chet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
682  $ ldu, tau, work, rwork, result( ntest ) )
683 *
684  CALL clacpy( ' ', n, n, v, ldu, a, lda )
685 *
686  ntest = ntest + 2
687  CALL cheevd_2stage( 'N', uplo, n, a, ldu, d3, work,
688  $ lwork, rwork, lrwedc, iwork, liwedc, iinfo )
689  IF( iinfo.NE.0 ) THEN
690  WRITE( nounit, fmt = 9999 )
691  $ 'CHEEVD_2STAGE(N,' // uplo //
692  $ ')', iinfo, n, jtype, ioldsd
693  info = abs( iinfo )
694  IF( iinfo.LT.0 ) THEN
695  RETURN
696  ELSE
697  result( ntest ) = ulpinv
698  GO TO 130
699  END IF
700  END IF
701 *
702 * Do test 3.
703 *
704  temp1 = zero
705  temp2 = zero
706  DO 120 j = 1, n
707  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
708  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
709  120 CONTINUE
710  result( ntest ) = temp2 / max( unfl,
711  $ ulp*max( temp1, temp2 ) )
712 *
713  130 CONTINUE
714  CALL clacpy( ' ', n, n, v, ldu, a, lda )
715 *
716  ntest = ntest + 1
717 *
718  IF( n.GT.0 ) THEN
719  temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
720  IF( il.NE.1 ) THEN
721  vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
722  $ ten*ulp*temp3, ten*rtunfl )
723  ELSE IF( n.GT.0 ) THEN
724  vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
725  $ ten*ulp*temp3, ten*rtunfl )
726  END IF
727  IF( iu.NE.n ) THEN
728  vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
729  $ ten*ulp*temp3, ten*rtunfl )
730  ELSE IF( n.GT.0 ) THEN
731  vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
732  $ ten*ulp*temp3, ten*rtunfl )
733  END IF
734  ELSE
735  temp3 = zero
736  vl = zero
737  vu = one
738  END IF
739 *
740  CALL cheevx( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
741  $ abstol, m, wa1, z, ldu, work, lwork, rwork,
742  $ iwork, iwork( 5*n+1 ), iinfo )
743  IF( iinfo.NE.0 ) THEN
744  WRITE( nounit, fmt = 9999 )'CHEEVX(V,A,' // uplo //
745  $ ')', iinfo, n, jtype, ioldsd
746  info = abs( iinfo )
747  IF( iinfo.LT.0 ) THEN
748  RETURN
749  ELSE
750  result( ntest ) = ulpinv
751  result( ntest+1 ) = ulpinv
752  result( ntest+2 ) = ulpinv
753  GO TO 150
754  END IF
755  END IF
756 *
757 * Do tests 4 and 5.
758 *
759  CALL clacpy( ' ', n, n, v, ldu, a, lda )
760 *
761  CALL chet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
762  $ ldu, tau, work, rwork, result( ntest ) )
763 *
764  ntest = ntest + 2
765  CALL cheevx_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
766  $ il, iu, abstol, m2, wa2, z, ldu,
767  $ work, lwork, rwork, iwork,
768  $ iwork( 5*n+1 ), iinfo )
769  IF( iinfo.NE.0 ) THEN
770  WRITE( nounit, fmt = 9999 )
771  $ 'CHEEVX_2STAGE(N,A,' // uplo //
772  $ ')', iinfo, n, jtype, ioldsd
773  info = abs( iinfo )
774  IF( iinfo.LT.0 ) THEN
775  RETURN
776  ELSE
777  result( ntest ) = ulpinv
778  GO TO 150
779  END IF
780  END IF
781 *
782 * Do test 6.
783 *
784  temp1 = zero
785  temp2 = zero
786  DO 140 j = 1, n
787  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
788  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
789  140 CONTINUE
790  result( ntest ) = temp2 / max( unfl,
791  $ ulp*max( temp1, temp2 ) )
792 *
793  150 CONTINUE
794  CALL clacpy( ' ', n, n, v, ldu, a, lda )
795 *
796  ntest = ntest + 1
797 *
798  CALL cheevx( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
799  $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
800  $ iwork, iwork( 5*n+1 ), iinfo )
801  IF( iinfo.NE.0 ) THEN
802  WRITE( nounit, fmt = 9999 )'CHEEVX(V,I,' // uplo //
803  $ ')', iinfo, n, jtype, ioldsd
804  info = abs( iinfo )
805  IF( iinfo.LT.0 ) THEN
806  RETURN
807  ELSE
808  result( ntest ) = ulpinv
809  GO TO 160
810  END IF
811  END IF
812 *
813 * Do tests 7 and 8.
814 *
815  CALL clacpy( ' ', n, n, v, ldu, a, lda )
816 *
817  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
818  $ v, ldu, tau, work, rwork, result( ntest ) )
819 *
820  ntest = ntest + 2
821 *
822  CALL cheevx_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
823  $ il, iu, abstol, m3, wa3, z, ldu,
824  $ work, lwork, rwork, iwork,
825  $ iwork( 5*n+1 ), iinfo )
826  IF( iinfo.NE.0 ) THEN
827  WRITE( nounit, fmt = 9999 )
828  $ 'CHEEVX_2STAGE(N,I,' // uplo //
829  $ ')', iinfo, n, jtype, ioldsd
830  info = abs( iinfo )
831  IF( iinfo.LT.0 ) THEN
832  RETURN
833  ELSE
834  result( ntest ) = ulpinv
835  GO TO 160
836  END IF
837  END IF
838 *
839 * Do test 9.
840 *
841  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
842  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
843  IF( n.GT.0 ) THEN
844  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
845  ELSE
846  temp3 = zero
847  END IF
848  result( ntest ) = ( temp1+temp2 ) /
849  $ max( unfl, temp3*ulp )
850 *
851  160 CONTINUE
852  CALL clacpy( ' ', n, n, v, ldu, a, lda )
853 *
854  ntest = ntest + 1
855 *
856  CALL cheevx( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
857  $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
858  $ iwork, iwork( 5*n+1 ), iinfo )
859  IF( iinfo.NE.0 ) THEN
860  WRITE( nounit, fmt = 9999 )'CHEEVX(V,V,' // uplo //
861  $ ')', iinfo, n, jtype, ioldsd
862  info = abs( iinfo )
863  IF( iinfo.LT.0 ) THEN
864  RETURN
865  ELSE
866  result( ntest ) = ulpinv
867  GO TO 170
868  END IF
869  END IF
870 *
871 * Do tests 10 and 11.
872 *
873  CALL clacpy( ' ', n, n, v, ldu, a, lda )
874 *
875  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
876  $ v, ldu, tau, work, rwork, result( ntest ) )
877 *
878  ntest = ntest + 2
879 *
880  CALL cheevx_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
881  $ il, iu, abstol, m3, wa3, z, ldu,
882  $ work, lwork, rwork, iwork,
883  $ iwork( 5*n+1 ), iinfo )
884  IF( iinfo.NE.0 ) THEN
885  WRITE( nounit, fmt = 9999 )
886  $ 'CHEEVX_2STAGE(N,V,' // uplo //
887  $ ')', iinfo, n, jtype, ioldsd
888  info = abs( iinfo )
889  IF( iinfo.LT.0 ) THEN
890  RETURN
891  ELSE
892  result( ntest ) = ulpinv
893  GO TO 170
894  END IF
895  END IF
896 *
897  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
898  result( ntest ) = ulpinv
899  GO TO 170
900  END IF
901 *
902 * Do test 12.
903 *
904  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
905  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
906  IF( n.GT.0 ) THEN
907  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
908  ELSE
909  temp3 = zero
910  END IF
911  result( ntest ) = ( temp1+temp2 ) /
912  $ max( unfl, temp3*ulp )
913 *
914  170 CONTINUE
915 *
916 * Call CHPEVD and CHPEVX.
917 *
918  CALL clacpy( ' ', n, n, v, ldu, a, lda )
919 *
920 * Load array WORK with the upper or lower triangular
921 * part of the matrix in packed form.
922 *
923  IF( iuplo.EQ.1 ) THEN
924  indx = 1
925  DO 190 j = 1, n
926  DO 180 i = 1, j
927  work( indx ) = a( i, j )
928  indx = indx + 1
929  180 CONTINUE
930  190 CONTINUE
931  ELSE
932  indx = 1
933  DO 210 j = 1, n
934  DO 200 i = j, n
935  work( indx ) = a( i, j )
936  indx = indx + 1
937  200 CONTINUE
938  210 CONTINUE
939  END IF
940 *
941  ntest = ntest + 1
942  indwrk = n*( n+1 ) / 2 + 1
943  CALL chpevd( 'V', uplo, n, work, d1, z, ldu,
944  $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
945  $ liwedc, iinfo )
946  IF( iinfo.NE.0 ) THEN
947  WRITE( nounit, fmt = 9999 )'CHPEVD(V,' // uplo //
948  $ ')', iinfo, n, jtype, ioldsd
949  info = abs( iinfo )
950  IF( iinfo.LT.0 ) THEN
951  RETURN
952  ELSE
953  result( ntest ) = ulpinv
954  result( ntest+1 ) = ulpinv
955  result( ntest+2 ) = ulpinv
956  GO TO 270
957  END IF
958  END IF
959 *
960 * Do tests 13 and 14.
961 *
962  CALL chet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
963  $ ldu, tau, work, rwork, result( ntest ) )
964 *
965  IF( iuplo.EQ.1 ) THEN
966  indx = 1
967  DO 230 j = 1, n
968  DO 220 i = 1, j
969  work( indx ) = a( i, j )
970  indx = indx + 1
971  220 CONTINUE
972  230 CONTINUE
973  ELSE
974  indx = 1
975  DO 250 j = 1, n
976  DO 240 i = j, n
977  work( indx ) = a( i, j )
978  indx = indx + 1
979  240 CONTINUE
980  250 CONTINUE
981  END IF
982 *
983  ntest = ntest + 2
984  indwrk = n*( n+1 ) / 2 + 1
985  CALL chpevd( 'N', uplo, n, work, d3, z, ldu,
986  $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
987  $ liwedc, iinfo )
988  IF( iinfo.NE.0 ) THEN
989  WRITE( nounit, fmt = 9999 )'CHPEVD(N,' // uplo //
990  $ ')', iinfo, n, jtype, ioldsd
991  info = abs( iinfo )
992  IF( iinfo.LT.0 ) THEN
993  RETURN
994  ELSE
995  result( ntest ) = ulpinv
996  GO TO 270
997  END IF
998  END IF
999 *
1000 * Do test 15.
1001 *
1002  temp1 = zero
1003  temp2 = zero
1004  DO 260 j = 1, n
1005  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1006  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1007  260 CONTINUE
1008  result( ntest ) = temp2 / max( unfl,
1009  $ ulp*max( temp1, temp2 ) )
1010 *
1011 * Load array WORK with the upper or lower triangular part
1012 * of the matrix in packed form.
1013 *
1014  270 CONTINUE
1015  IF( iuplo.EQ.1 ) THEN
1016  indx = 1
1017  DO 290 j = 1, n
1018  DO 280 i = 1, j
1019  work( indx ) = a( i, j )
1020  indx = indx + 1
1021  280 CONTINUE
1022  290 CONTINUE
1023  ELSE
1024  indx = 1
1025  DO 310 j = 1, n
1026  DO 300 i = j, n
1027  work( indx ) = a( i, j )
1028  indx = indx + 1
1029  300 CONTINUE
1030  310 CONTINUE
1031  END IF
1032 *
1033  ntest = ntest + 1
1034 *
1035  IF( n.GT.0 ) THEN
1036  temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1037  IF( il.NE.1 ) THEN
1038  vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1039  $ ten*ulp*temp3, ten*rtunfl )
1040  ELSE IF( n.GT.0 ) THEN
1041  vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1042  $ ten*ulp*temp3, ten*rtunfl )
1043  END IF
1044  IF( iu.NE.n ) THEN
1045  vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1046  $ ten*ulp*temp3, ten*rtunfl )
1047  ELSE IF( n.GT.0 ) THEN
1048  vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1049  $ ten*ulp*temp3, ten*rtunfl )
1050  END IF
1051  ELSE
1052  temp3 = zero
1053  vl = zero
1054  vu = one
1055  END IF
1056 *
1057  CALL chpevx( 'V', 'A', uplo, n, work, vl, vu, il, iu,
1058  $ abstol, m, wa1, z, ldu, v, rwork, iwork,
1059  $ iwork( 5*n+1 ), iinfo )
1060  IF( iinfo.NE.0 ) THEN
1061  WRITE( nounit, fmt = 9999 )'CHPEVX(V,A,' // uplo //
1062  $ ')', iinfo, n, jtype, ioldsd
1063  info = abs( iinfo )
1064  IF( iinfo.LT.0 ) THEN
1065  RETURN
1066  ELSE
1067  result( ntest ) = ulpinv
1068  result( ntest+1 ) = ulpinv
1069  result( ntest+2 ) = ulpinv
1070  GO TO 370
1071  END IF
1072  END IF
1073 *
1074 * Do tests 16 and 17.
1075 *
1076  CALL chet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1077  $ ldu, tau, work, rwork, result( ntest ) )
1078 *
1079  ntest = ntest + 2
1080 *
1081  IF( iuplo.EQ.1 ) THEN
1082  indx = 1
1083  DO 330 j = 1, n
1084  DO 320 i = 1, j
1085  work( indx ) = a( i, j )
1086  indx = indx + 1
1087  320 CONTINUE
1088  330 CONTINUE
1089  ELSE
1090  indx = 1
1091  DO 350 j = 1, n
1092  DO 340 i = j, n
1093  work( indx ) = a( i, j )
1094  indx = indx + 1
1095  340 CONTINUE
1096  350 CONTINUE
1097  END IF
1098 *
1099  CALL chpevx( 'N', 'A', uplo, n, work, vl, vu, il, iu,
1100  $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1101  $ iwork( 5*n+1 ), iinfo )
1102  IF( iinfo.NE.0 ) THEN
1103  WRITE( nounit, fmt = 9999 )'CHPEVX(N,A,' // uplo //
1104  $ ')', iinfo, n, jtype, ioldsd
1105  info = abs( iinfo )
1106  IF( iinfo.LT.0 ) THEN
1107  RETURN
1108  ELSE
1109  result( ntest ) = ulpinv
1110  GO TO 370
1111  END IF
1112  END IF
1113 *
1114 * Do test 18.
1115 *
1116  temp1 = zero
1117  temp2 = zero
1118  DO 360 j = 1, n
1119  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1120  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1121  360 CONTINUE
1122  result( ntest ) = temp2 / max( unfl,
1123  $ ulp*max( temp1, temp2 ) )
1124 *
1125  370 CONTINUE
1126  ntest = ntest + 1
1127  IF( iuplo.EQ.1 ) THEN
1128  indx = 1
1129  DO 390 j = 1, n
1130  DO 380 i = 1, j
1131  work( indx ) = a( i, j )
1132  indx = indx + 1
1133  380 CONTINUE
1134  390 CONTINUE
1135  ELSE
1136  indx = 1
1137  DO 410 j = 1, n
1138  DO 400 i = j, n
1139  work( indx ) = a( i, j )
1140  indx = indx + 1
1141  400 CONTINUE
1142  410 CONTINUE
1143  END IF
1144 *
1145  CALL chpevx( 'V', 'I', uplo, n, work, vl, vu, il, iu,
1146  $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1147  $ iwork( 5*n+1 ), iinfo )
1148  IF( iinfo.NE.0 ) THEN
1149  WRITE( nounit, fmt = 9999 )'CHPEVX(V,I,' // uplo //
1150  $ ')', iinfo, n, jtype, ioldsd
1151  info = abs( iinfo )
1152  IF( iinfo.LT.0 ) THEN
1153  RETURN
1154  ELSE
1155  result( ntest ) = ulpinv
1156  result( ntest+1 ) = ulpinv
1157  result( ntest+2 ) = ulpinv
1158  GO TO 460
1159  END IF
1160  END IF
1161 *
1162 * Do tests 19 and 20.
1163 *
1164  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1165  $ v, ldu, tau, work, rwork, result( ntest ) )
1166 *
1167  ntest = ntest + 2
1168 *
1169  IF( iuplo.EQ.1 ) THEN
1170  indx = 1
1171  DO 430 j = 1, n
1172  DO 420 i = 1, j
1173  work( indx ) = a( i, j )
1174  indx = indx + 1
1175  420 CONTINUE
1176  430 CONTINUE
1177  ELSE
1178  indx = 1
1179  DO 450 j = 1, n
1180  DO 440 i = j, n
1181  work( indx ) = a( i, j )
1182  indx = indx + 1
1183  440 CONTINUE
1184  450 CONTINUE
1185  END IF
1186 *
1187  CALL chpevx( 'N', 'I', uplo, n, work, vl, vu, il, iu,
1188  $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1189  $ iwork( 5*n+1 ), iinfo )
1190  IF( iinfo.NE.0 ) THEN
1191  WRITE( nounit, fmt = 9999 )'CHPEVX(N,I,' // uplo //
1192  $ ')', iinfo, n, jtype, ioldsd
1193  info = abs( iinfo )
1194  IF( iinfo.LT.0 ) THEN
1195  RETURN
1196  ELSE
1197  result( ntest ) = ulpinv
1198  GO TO 460
1199  END IF
1200  END IF
1201 *
1202 * Do test 21.
1203 *
1204  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1205  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1206  IF( n.GT.0 ) THEN
1207  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1208  ELSE
1209  temp3 = zero
1210  END IF
1211  result( ntest ) = ( temp1+temp2 ) /
1212  $ max( unfl, temp3*ulp )
1213 *
1214  460 CONTINUE
1215  ntest = ntest + 1
1216  IF( iuplo.EQ.1 ) THEN
1217  indx = 1
1218  DO 480 j = 1, n
1219  DO 470 i = 1, j
1220  work( indx ) = a( i, j )
1221  indx = indx + 1
1222  470 CONTINUE
1223  480 CONTINUE
1224  ELSE
1225  indx = 1
1226  DO 500 j = 1, n
1227  DO 490 i = j, n
1228  work( indx ) = a( i, j )
1229  indx = indx + 1
1230  490 CONTINUE
1231  500 CONTINUE
1232  END IF
1233 *
1234  CALL chpevx( 'V', 'V', uplo, n, work, vl, vu, il, iu,
1235  $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1236  $ iwork( 5*n+1 ), iinfo )
1237  IF( iinfo.NE.0 ) THEN
1238  WRITE( nounit, fmt = 9999 )'CHPEVX(V,V,' // uplo //
1239  $ ')', iinfo, n, jtype, ioldsd
1240  info = abs( iinfo )
1241  IF( iinfo.LT.0 ) THEN
1242  RETURN
1243  ELSE
1244  result( ntest ) = ulpinv
1245  result( ntest+1 ) = ulpinv
1246  result( ntest+2 ) = ulpinv
1247  GO TO 550
1248  END IF
1249  END IF
1250 *
1251 * Do tests 22 and 23.
1252 *
1253  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1254  $ v, ldu, tau, work, rwork, result( ntest ) )
1255 *
1256  ntest = ntest + 2
1257 *
1258  IF( iuplo.EQ.1 ) THEN
1259  indx = 1
1260  DO 520 j = 1, n
1261  DO 510 i = 1, j
1262  work( indx ) = a( i, j )
1263  indx = indx + 1
1264  510 CONTINUE
1265  520 CONTINUE
1266  ELSE
1267  indx = 1
1268  DO 540 j = 1, n
1269  DO 530 i = j, n
1270  work( indx ) = a( i, j )
1271  indx = indx + 1
1272  530 CONTINUE
1273  540 CONTINUE
1274  END IF
1275 *
1276  CALL chpevx( 'N', 'V', uplo, n, work, vl, vu, il, iu,
1277  $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1278  $ iwork( 5*n+1 ), iinfo )
1279  IF( iinfo.NE.0 ) THEN
1280  WRITE( nounit, fmt = 9999 )'CHPEVX(N,V,' // uplo //
1281  $ ')', iinfo, n, jtype, ioldsd
1282  info = abs( iinfo )
1283  IF( iinfo.LT.0 ) THEN
1284  RETURN
1285  ELSE
1286  result( ntest ) = ulpinv
1287  GO TO 550
1288  END IF
1289  END IF
1290 *
1291  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1292  result( ntest ) = ulpinv
1293  GO TO 550
1294  END IF
1295 *
1296 * Do test 24.
1297 *
1298  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1299  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1300  IF( n.GT.0 ) THEN
1301  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1302  ELSE
1303  temp3 = zero
1304  END IF
1305  result( ntest ) = ( temp1+temp2 ) /
1306  $ max( unfl, temp3*ulp )
1307 *
1308  550 CONTINUE
1309 *
1310 * Call CHBEVD and CHBEVX.
1311 *
1312  IF( jtype.LE.7 ) THEN
1313  kd = 0
1314  ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1315  kd = max( n-1, 0 )
1316  ELSE
1317  kd = ihbw
1318  END IF
1319 *
1320 * Load array V with the upper or lower triangular part
1321 * of the matrix in band form.
1322 *
1323  IF( iuplo.EQ.1 ) THEN
1324  DO 570 j = 1, n
1325  DO 560 i = max( 1, j-kd ), j
1326  v( kd+1+i-j, j ) = a( i, j )
1327  560 CONTINUE
1328  570 CONTINUE
1329  ELSE
1330  DO 590 j = 1, n
1331  DO 580 i = j, min( n, j+kd )
1332  v( 1+i-j, j ) = a( i, j )
1333  580 CONTINUE
1334  590 CONTINUE
1335  END IF
1336 *
1337  ntest = ntest + 1
1338  CALL chbevd( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1339  $ lwedc, rwork, lrwedc, iwork, liwedc, iinfo )
1340  IF( iinfo.NE.0 ) THEN
1341  WRITE( nounit, fmt = 9998 )'CHBEVD(V,' // uplo //
1342  $ ')', iinfo, n, kd, jtype, ioldsd
1343  info = abs( iinfo )
1344  IF( iinfo.LT.0 ) THEN
1345  RETURN
1346  ELSE
1347  result( ntest ) = ulpinv
1348  result( ntest+1 ) = ulpinv
1349  result( ntest+2 ) = ulpinv
1350  GO TO 650
1351  END IF
1352  END IF
1353 *
1354 * Do tests 25 and 26.
1355 *
1356  CALL chet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1357  $ ldu, tau, work, rwork, result( ntest ) )
1358 *
1359  IF( iuplo.EQ.1 ) THEN
1360  DO 610 j = 1, n
1361  DO 600 i = max( 1, j-kd ), j
1362  v( kd+1+i-j, j ) = a( i, j )
1363  600 CONTINUE
1364  610 CONTINUE
1365  ELSE
1366  DO 630 j = 1, n
1367  DO 620 i = j, min( n, j+kd )
1368  v( 1+i-j, j ) = a( i, j )
1369  620 CONTINUE
1370  630 CONTINUE
1371  END IF
1372 *
1373  ntest = ntest + 2
1374  CALL chbevd_2stage( 'N', uplo, n, kd, v, ldu, d3,
1375  $ z, ldu, work, lwork, rwork,
1376  $ lrwedc, iwork, liwedc, iinfo )
1377  IF( iinfo.NE.0 ) THEN
1378  WRITE( nounit, fmt = 9998 )
1379  $ 'CHBEVD_2STAGE(N,' // uplo //
1380  $ ')', iinfo, n, kd, jtype, ioldsd
1381  info = abs( iinfo )
1382  IF( iinfo.LT.0 ) THEN
1383  RETURN
1384  ELSE
1385  result( ntest ) = ulpinv
1386  GO TO 650
1387  END IF
1388  END IF
1389 *
1390 * Do test 27.
1391 *
1392  temp1 = zero
1393  temp2 = zero
1394  DO 640 j = 1, n
1395  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1396  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1397  640 CONTINUE
1398  result( ntest ) = temp2 / max( unfl,
1399  $ ulp*max( temp1, temp2 ) )
1400 *
1401 * Load array V with the upper or lower triangular part
1402 * of the matrix in band form.
1403 *
1404  650 CONTINUE
1405  IF( iuplo.EQ.1 ) THEN
1406  DO 670 j = 1, n
1407  DO 660 i = max( 1, j-kd ), j
1408  v( kd+1+i-j, j ) = a( i, j )
1409  660 CONTINUE
1410  670 CONTINUE
1411  ELSE
1412  DO 690 j = 1, n
1413  DO 680 i = j, min( n, j+kd )
1414  v( 1+i-j, j ) = a( i, j )
1415  680 CONTINUE
1416  690 CONTINUE
1417  END IF
1418 *
1419  ntest = ntest + 1
1420  CALL chbevx( 'V', 'A', uplo, n, kd, v, ldu, u, ldu, vl,
1421  $ vu, il, iu, abstol, m, wa1, z, ldu, work,
1422  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1423  IF( iinfo.NE.0 ) THEN
1424  WRITE( nounit, fmt = 9999 )'CHBEVX(V,A,' // uplo //
1425  $ ')', iinfo, n, kd, jtype, ioldsd
1426  info = abs( iinfo )
1427  IF( iinfo.LT.0 ) THEN
1428  RETURN
1429  ELSE
1430  result( ntest ) = ulpinv
1431  result( ntest+1 ) = ulpinv
1432  result( ntest+2 ) = ulpinv
1433  GO TO 750
1434  END IF
1435  END IF
1436 *
1437 * Do tests 28 and 29.
1438 *
1439  CALL chet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1440  $ ldu, tau, work, rwork, result( ntest ) )
1441 *
1442  ntest = ntest + 2
1443 *
1444  IF( iuplo.EQ.1 ) THEN
1445  DO 710 j = 1, n
1446  DO 700 i = max( 1, j-kd ), j
1447  v( kd+1+i-j, j ) = a( i, j )
1448  700 CONTINUE
1449  710 CONTINUE
1450  ELSE
1451  DO 730 j = 1, n
1452  DO 720 i = j, min( n, j+kd )
1453  v( 1+i-j, j ) = a( i, j )
1454  720 CONTINUE
1455  730 CONTINUE
1456  END IF
1457 *
1458  CALL chbevx_2stage( 'N', 'A', uplo, n, kd, v, ldu,
1459  $ u, ldu, vl, vu, il, iu, abstol,
1460  $ m2, wa2, z, ldu, work, lwork,
1461  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1462  IF( iinfo.NE.0 ) THEN
1463  WRITE( nounit, fmt = 9998 )
1464  $ 'CHBEVX_2STAGE(N,A,' // uplo //
1465  $ ')', iinfo, n, kd, jtype, ioldsd
1466  info = abs( iinfo )
1467  IF( iinfo.LT.0 ) THEN
1468  RETURN
1469  ELSE
1470  result( ntest ) = ulpinv
1471  GO TO 750
1472  END IF
1473  END IF
1474 *
1475 * Do test 30.
1476 *
1477  temp1 = zero
1478  temp2 = zero
1479  DO 740 j = 1, n
1480  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1481  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1482  740 CONTINUE
1483  result( ntest ) = temp2 / max( unfl,
1484  $ ulp*max( temp1, temp2 ) )
1485 *
1486 * Load array V with the upper or lower triangular part
1487 * of the matrix in band form.
1488 *
1489  750 CONTINUE
1490  ntest = ntest + 1
1491  IF( iuplo.EQ.1 ) THEN
1492  DO 770 j = 1, n
1493  DO 760 i = max( 1, j-kd ), j
1494  v( kd+1+i-j, j ) = a( i, j )
1495  760 CONTINUE
1496  770 CONTINUE
1497  ELSE
1498  DO 790 j = 1, n
1499  DO 780 i = j, min( n, j+kd )
1500  v( 1+i-j, j ) = a( i, j )
1501  780 CONTINUE
1502  790 CONTINUE
1503  END IF
1504 *
1505  CALL chbevx( 'V', 'I', uplo, n, kd, v, ldu, u, ldu, vl,
1506  $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1507  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1508  IF( iinfo.NE.0 ) THEN
1509  WRITE( nounit, fmt = 9998 )'CHBEVX(V,I,' // uplo //
1510  $ ')', iinfo, n, kd, jtype, ioldsd
1511  info = abs( iinfo )
1512  IF( iinfo.LT.0 ) THEN
1513  RETURN
1514  ELSE
1515  result( ntest ) = ulpinv
1516  result( ntest+1 ) = ulpinv
1517  result( ntest+2 ) = ulpinv
1518  GO TO 840
1519  END IF
1520  END IF
1521 *
1522 * Do tests 31 and 32.
1523 *
1524  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1525  $ v, ldu, tau, work, rwork, result( ntest ) )
1526 *
1527  ntest = ntest + 2
1528 *
1529  IF( iuplo.EQ.1 ) THEN
1530  DO 810 j = 1, n
1531  DO 800 i = max( 1, j-kd ), j
1532  v( kd+1+i-j, j ) = a( i, j )
1533  800 CONTINUE
1534  810 CONTINUE
1535  ELSE
1536  DO 830 j = 1, n
1537  DO 820 i = j, min( n, j+kd )
1538  v( 1+i-j, j ) = a( i, j )
1539  820 CONTINUE
1540  830 CONTINUE
1541  END IF
1542  CALL chbevx_2stage( 'N', 'I', uplo, n, kd, v, ldu,
1543  $ u, ldu, vl, vu, il, iu, abstol,
1544  $ m3, wa3, z, ldu, work, lwork,
1545  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1546  IF( iinfo.NE.0 ) THEN
1547  WRITE( nounit, fmt = 9998 )
1548  $ 'CHBEVX_2STAGE(N,I,' // uplo //
1549  $ ')', iinfo, n, kd, jtype, ioldsd
1550  info = abs( iinfo )
1551  IF( iinfo.LT.0 ) THEN
1552  RETURN
1553  ELSE
1554  result( ntest ) = ulpinv
1555  GO TO 840
1556  END IF
1557  END IF
1558 *
1559 * Do test 33.
1560 *
1561  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1562  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1563  IF( n.GT.0 ) THEN
1564  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1565  ELSE
1566  temp3 = zero
1567  END IF
1568  result( ntest ) = ( temp1+temp2 ) /
1569  $ max( unfl, temp3*ulp )
1570 *
1571 * Load array V with the upper or lower triangular part
1572 * of the matrix in band form.
1573 *
1574  840 CONTINUE
1575  ntest = ntest + 1
1576  IF( iuplo.EQ.1 ) THEN
1577  DO 860 j = 1, n
1578  DO 850 i = max( 1, j-kd ), j
1579  v( kd+1+i-j, j ) = a( i, j )
1580  850 CONTINUE
1581  860 CONTINUE
1582  ELSE
1583  DO 880 j = 1, n
1584  DO 870 i = j, min( n, j+kd )
1585  v( 1+i-j, j ) = a( i, j )
1586  870 CONTINUE
1587  880 CONTINUE
1588  END IF
1589  CALL chbevx( 'V', 'V', uplo, n, kd, v, ldu, u, ldu, vl,
1590  $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1591  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1592  IF( iinfo.NE.0 ) THEN
1593  WRITE( nounit, fmt = 9998 )'CHBEVX(V,V,' // uplo //
1594  $ ')', iinfo, n, kd, jtype, ioldsd
1595  info = abs( iinfo )
1596  IF( iinfo.LT.0 ) THEN
1597  RETURN
1598  ELSE
1599  result( ntest ) = ulpinv
1600  result( ntest+1 ) = ulpinv
1601  result( ntest+2 ) = ulpinv
1602  GO TO 930
1603  END IF
1604  END IF
1605 *
1606 * Do tests 34 and 35.
1607 *
1608  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1609  $ v, ldu, tau, work, rwork, result( ntest ) )
1610 *
1611  ntest = ntest + 2
1612 *
1613  IF( iuplo.EQ.1 ) THEN
1614  DO 900 j = 1, n
1615  DO 890 i = max( 1, j-kd ), j
1616  v( kd+1+i-j, j ) = a( i, j )
1617  890 CONTINUE
1618  900 CONTINUE
1619  ELSE
1620  DO 920 j = 1, n
1621  DO 910 i = j, min( n, j+kd )
1622  v( 1+i-j, j ) = a( i, j )
1623  910 CONTINUE
1624  920 CONTINUE
1625  END IF
1626  CALL chbevx_2stage( 'N', 'V', uplo, n, kd, v, ldu,
1627  $ u, ldu, vl, vu, il, iu, abstol,
1628  $ m3, wa3, z, ldu, work, lwork,
1629  $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1630  IF( iinfo.NE.0 ) THEN
1631  WRITE( nounit, fmt = 9998 )
1632  $ 'CHBEVX_2STAGE(N,V,' // uplo //
1633  $ ')', iinfo, n, kd, jtype, ioldsd
1634  info = abs( iinfo )
1635  IF( iinfo.LT.0 ) THEN
1636  RETURN
1637  ELSE
1638  result( ntest ) = ulpinv
1639  GO TO 930
1640  END IF
1641  END IF
1642 *
1643  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1644  result( ntest ) = ulpinv
1645  GO TO 930
1646  END IF
1647 *
1648 * Do test 36.
1649 *
1650  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1651  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1652  IF( n.GT.0 ) THEN
1653  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1654  ELSE
1655  temp3 = zero
1656  END IF
1657  result( ntest ) = ( temp1+temp2 ) /
1658  $ max( unfl, temp3*ulp )
1659 *
1660  930 CONTINUE
1661 *
1662 * Call CHEEV
1663 *
1664  CALL clacpy( ' ', n, n, a, lda, v, ldu )
1665 *
1666  ntest = ntest + 1
1667  CALL cheev( 'V', uplo, n, a, ldu, d1, work, lwork, rwork,
1668  $ iinfo )
1669  IF( iinfo.NE.0 ) THEN
1670  WRITE( nounit, fmt = 9999 )'CHEEV(V,' // uplo // ')',
1671  $ iinfo, n, jtype, ioldsd
1672  info = abs( iinfo )
1673  IF( iinfo.LT.0 ) THEN
1674  RETURN
1675  ELSE
1676  result( ntest ) = ulpinv
1677  result( ntest+1 ) = ulpinv
1678  result( ntest+2 ) = ulpinv
1679  GO TO 950
1680  END IF
1681  END IF
1682 *
1683 * Do tests 37 and 38
1684 *
1685  CALL chet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
1686  $ ldu, tau, work, rwork, result( ntest ) )
1687 *
1688  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1689 *
1690  ntest = ntest + 2
1691  CALL cheev_2stage( 'N', uplo, n, a, ldu, d3,
1692  $ work, lwork, rwork, iinfo )
1693  IF( iinfo.NE.0 ) THEN
1694  WRITE( nounit, fmt = 9999 )
1695  $ 'CHEEV_2STAGE(N,' // uplo // ')',
1696  $ iinfo, n, jtype, ioldsd
1697  info = abs( iinfo )
1698  IF( iinfo.LT.0 ) THEN
1699  RETURN
1700  ELSE
1701  result( ntest ) = ulpinv
1702  GO TO 950
1703  END IF
1704  END IF
1705 *
1706 * Do test 39
1707 *
1708  temp1 = zero
1709  temp2 = zero
1710  DO 940 j = 1, n
1711  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1712  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1713  940 CONTINUE
1714  result( ntest ) = temp2 / max( unfl,
1715  $ ulp*max( temp1, temp2 ) )
1716 *
1717  950 CONTINUE
1718 *
1719  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1720 *
1721 * Call CHPEV
1722 *
1723 * Load array WORK with the upper or lower triangular
1724 * part of the matrix in packed form.
1725 *
1726  IF( iuplo.EQ.1 ) THEN
1727  indx = 1
1728  DO 970 j = 1, n
1729  DO 960 i = 1, j
1730  work( indx ) = a( i, j )
1731  indx = indx + 1
1732  960 CONTINUE
1733  970 CONTINUE
1734  ELSE
1735  indx = 1
1736  DO 990 j = 1, n
1737  DO 980 i = j, n
1738  work( indx ) = a( i, j )
1739  indx = indx + 1
1740  980 CONTINUE
1741  990 CONTINUE
1742  END IF
1743 *
1744  ntest = ntest + 1
1745  indwrk = n*( n+1 ) / 2 + 1
1746  CALL chpev( 'V', uplo, n, work, d1, z, ldu,
1747  $ work( indwrk ), rwork, iinfo )
1748  IF( iinfo.NE.0 ) THEN
1749  WRITE( nounit, fmt = 9999 )'CHPEV(V,' // uplo // ')',
1750  $ iinfo, n, jtype, ioldsd
1751  info = abs( iinfo )
1752  IF( iinfo.LT.0 ) THEN
1753  RETURN
1754  ELSE
1755  result( ntest ) = ulpinv
1756  result( ntest+1 ) = ulpinv
1757  result( ntest+2 ) = ulpinv
1758  GO TO 1050
1759  END IF
1760  END IF
1761 *
1762 * Do tests 40 and 41.
1763 *
1764  CALL chet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1765  $ ldu, tau, work, rwork, result( ntest ) )
1766 *
1767  IF( iuplo.EQ.1 ) THEN
1768  indx = 1
1769  DO 1010 j = 1, n
1770  DO 1000 i = 1, j
1771  work( indx ) = a( i, j )
1772  indx = indx + 1
1773  1000 CONTINUE
1774  1010 CONTINUE
1775  ELSE
1776  indx = 1
1777  DO 1030 j = 1, n
1778  DO 1020 i = j, n
1779  work( indx ) = a( i, j )
1780  indx = indx + 1
1781  1020 CONTINUE
1782  1030 CONTINUE
1783  END IF
1784 *
1785  ntest = ntest + 2
1786  indwrk = n*( n+1 ) / 2 + 1
1787  CALL chpev( 'N', uplo, n, work, d3, z, ldu,
1788  $ work( indwrk ), rwork, iinfo )
1789  IF( iinfo.NE.0 ) THEN
1790  WRITE( nounit, fmt = 9999 )'CHPEV(N,' // uplo // ')',
1791  $ iinfo, n, jtype, ioldsd
1792  info = abs( iinfo )
1793  IF( iinfo.LT.0 ) THEN
1794  RETURN
1795  ELSE
1796  result( ntest ) = ulpinv
1797  GO TO 1050
1798  END IF
1799  END IF
1800 *
1801 * Do test 42
1802 *
1803  temp1 = zero
1804  temp2 = zero
1805  DO 1040 j = 1, n
1806  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1807  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1808  1040 CONTINUE
1809  result( ntest ) = temp2 / max( unfl,
1810  $ ulp*max( temp1, temp2 ) )
1811 *
1812  1050 CONTINUE
1813 *
1814 * Call CHBEV
1815 *
1816  IF( jtype.LE.7 ) THEN
1817  kd = 0
1818  ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1819  kd = max( n-1, 0 )
1820  ELSE
1821  kd = ihbw
1822  END IF
1823 *
1824 * Load array V with the upper or lower triangular part
1825 * of the matrix in band form.
1826 *
1827  IF( iuplo.EQ.1 ) THEN
1828  DO 1070 j = 1, n
1829  DO 1060 i = max( 1, j-kd ), j
1830  v( kd+1+i-j, j ) = a( i, j )
1831  1060 CONTINUE
1832  1070 CONTINUE
1833  ELSE
1834  DO 1090 j = 1, n
1835  DO 1080 i = j, min( n, j+kd )
1836  v( 1+i-j, j ) = a( i, j )
1837  1080 CONTINUE
1838  1090 CONTINUE
1839  END IF
1840 *
1841  ntest = ntest + 1
1842  CALL chbev( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1843  $ rwork, iinfo )
1844  IF( iinfo.NE.0 ) THEN
1845  WRITE( nounit, fmt = 9998 )'CHBEV(V,' // uplo // ')',
1846  $ iinfo, n, kd, jtype, ioldsd
1847  info = abs( iinfo )
1848  IF( iinfo.LT.0 ) THEN
1849  RETURN
1850  ELSE
1851  result( ntest ) = ulpinv
1852  result( ntest+1 ) = ulpinv
1853  result( ntest+2 ) = ulpinv
1854  GO TO 1140
1855  END IF
1856  END IF
1857 *
1858 * Do tests 43 and 44.
1859 *
1860  CALL chet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1861  $ ldu, tau, work, rwork, result( ntest ) )
1862 *
1863  IF( iuplo.EQ.1 ) THEN
1864  DO 1110 j = 1, n
1865  DO 1100 i = max( 1, j-kd ), j
1866  v( kd+1+i-j, j ) = a( i, j )
1867  1100 CONTINUE
1868  1110 CONTINUE
1869  ELSE
1870  DO 1130 j = 1, n
1871  DO 1120 i = j, min( n, j+kd )
1872  v( 1+i-j, j ) = a( i, j )
1873  1120 CONTINUE
1874  1130 CONTINUE
1875  END IF
1876 *
1877  ntest = ntest + 2
1878  CALL chbev_2stage( 'N', uplo, n, kd, v, ldu, d3, z, ldu,
1879  $ work, lwork, rwork, iinfo )
1880  IF( iinfo.NE.0 ) THEN
1881  WRITE( nounit, fmt = 9998 )
1882  $ 'CHBEV_2STAGE(N,' // uplo // ')',
1883  $ iinfo, n, kd, jtype, ioldsd
1884  info = abs( iinfo )
1885  IF( iinfo.LT.0 ) THEN
1886  RETURN
1887  ELSE
1888  result( ntest ) = ulpinv
1889  GO TO 1140
1890  END IF
1891  END IF
1892 *
1893  1140 CONTINUE
1894 *
1895 * Do test 45.
1896 *
1897  temp1 = zero
1898  temp2 = zero
1899  DO 1150 j = 1, n
1900  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1901  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1902  1150 CONTINUE
1903  result( ntest ) = temp2 / max( unfl,
1904  $ ulp*max( temp1, temp2 ) )
1905 *
1906  CALL clacpy( ' ', n, n, a, lda, v, ldu )
1907  ntest = ntest + 1
1908  CALL cheevr( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
1909  $ abstol, m, wa1, z, ldu, iwork, work, lwork,
1910  $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1911  $ iinfo )
1912  IF( iinfo.NE.0 ) THEN
1913  WRITE( nounit, fmt = 9999 )'CHEEVR(V,A,' // uplo //
1914  $ ')', iinfo, n, jtype, ioldsd
1915  info = abs( iinfo )
1916  IF( iinfo.LT.0 ) THEN
1917  RETURN
1918  ELSE
1919  result( ntest ) = ulpinv
1920  result( ntest+1 ) = ulpinv
1921  result( ntest+2 ) = ulpinv
1922  GO TO 1170
1923  END IF
1924  END IF
1925 *
1926 * Do tests 45 and 46 (or ... )
1927 *
1928  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1929 *
1930  CALL chet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1931  $ ldu, tau, work, rwork, result( ntest ) )
1932 *
1933  ntest = ntest + 2
1934  CALL cheevr_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
1935  $ il, iu, abstol, m2, wa2, z, ldu,
1936  $ iwork, work, lwork, rwork, lrwork,
1937  $ iwork( 2*n+1 ), liwork-2*n, iinfo )
1938  IF( iinfo.NE.0 ) THEN
1939  WRITE( nounit, fmt = 9999 )
1940  $ 'CHEEVR_2STAGE(N,A,' // uplo //
1941  $ ')', iinfo, n, jtype, ioldsd
1942  info = abs( iinfo )
1943  IF( iinfo.LT.0 ) THEN
1944  RETURN
1945  ELSE
1946  result( ntest ) = ulpinv
1947  GO TO 1170
1948  END IF
1949  END IF
1950 *
1951 * Do test 47 (or ... )
1952 *
1953  temp1 = zero
1954  temp2 = zero
1955  DO 1160 j = 1, n
1956  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1957  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1958  1160 CONTINUE
1959  result( ntest ) = temp2 / max( unfl,
1960  $ ulp*max( temp1, temp2 ) )
1961 *
1962  1170 CONTINUE
1963 *
1964  ntest = ntest + 1
1965  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1966  CALL cheevr( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
1967  $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
1968  $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1969  $ iinfo )
1970  IF( iinfo.NE.0 ) THEN
1971  WRITE( nounit, fmt = 9999 )'CHEEVR(V,I,' // uplo //
1972  $ ')', iinfo, n, jtype, ioldsd
1973  info = abs( iinfo )
1974  IF( iinfo.LT.0 ) THEN
1975  RETURN
1976  ELSE
1977  result( ntest ) = ulpinv
1978  result( ntest+1 ) = ulpinv
1979  result( ntest+2 ) = ulpinv
1980  GO TO 1180
1981  END IF
1982  END IF
1983 *
1984 * Do tests 48 and 49 (or +??)
1985 *
1986  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1987 *
1988  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1989  $ v, ldu, tau, work, rwork, result( ntest ) )
1990 *
1991  ntest = ntest + 2
1992  CALL clacpy( ' ', n, n, v, ldu, a, lda )
1993  CALL cheevr_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
1994  $ il, iu, abstol, m3, wa3, z, ldu,
1995  $ iwork, work, lwork, rwork, lrwork,
1996  $ iwork( 2*n+1 ), liwork-2*n, iinfo )
1997  IF( iinfo.NE.0 ) THEN
1998  WRITE( nounit, fmt = 9999 )
1999  $ 'CHEEVR_2STAGE(N,I,' // uplo //
2000  $ ')', iinfo, n, jtype, ioldsd
2001  info = abs( iinfo )
2002  IF( iinfo.LT.0 ) THEN
2003  RETURN
2004  ELSE
2005  result( ntest ) = ulpinv
2006  GO TO 1180
2007  END IF
2008  END IF
2009 *
2010 * Do test 50 (or +??)
2011 *
2012  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2013  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2014  result( ntest ) = ( temp1+temp2 ) /
2015  $ max( unfl, ulp*temp3 )
2016  1180 CONTINUE
2017 *
2018  ntest = ntest + 1
2019  CALL clacpy( ' ', n, n, v, ldu, a, lda )
2020  CALL cheevr( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
2021  $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2022  $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
2023  $ iinfo )
2024  IF( iinfo.NE.0 ) THEN
2025  WRITE( nounit, fmt = 9999 )'CHEEVR(V,V,' // uplo //
2026  $ ')', iinfo, n, jtype, ioldsd
2027  info = abs( iinfo )
2028  IF( iinfo.LT.0 ) THEN
2029  RETURN
2030  ELSE
2031  result( ntest ) = ulpinv
2032  result( ntest+1 ) = ulpinv
2033  result( ntest+2 ) = ulpinv
2034  GO TO 1190
2035  END IF
2036  END IF
2037 *
2038 * Do tests 51 and 52 (or +??)
2039 *
2040  CALL clacpy( ' ', n, n, v, ldu, a, lda )
2041 *
2042  CALL chet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2043  $ v, ldu, tau, work, rwork, result( ntest ) )
2044 *
2045  ntest = ntest + 2
2046  CALL clacpy( ' ', n, n, v, ldu, a, lda )
2047  CALL cheevr_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
2048  $ il, iu, abstol, m3, wa3, z, ldu,
2049  $ iwork, work, lwork, rwork, lrwork,
2050  $ iwork( 2*n+1 ), liwork-2*n, iinfo )
2051  IF( iinfo.NE.0 ) THEN
2052  WRITE( nounit, fmt = 9999 )
2053  $ 'CHEEVR_2STAGE(N,V,' // uplo //
2054  $ ')', iinfo, n, jtype, ioldsd
2055  info = abs( iinfo )
2056  IF( iinfo.LT.0 ) THEN
2057  RETURN
2058  ELSE
2059  result( ntest ) = ulpinv
2060  GO TO 1190
2061  END IF
2062  END IF
2063 *
2064  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2065  result( ntest ) = ulpinv
2066  GO TO 1190
2067  END IF
2068 *
2069 * Do test 52 (or +??)
2070 *
2071  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2072  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2073  IF( n.GT.0 ) THEN
2074  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2075  ELSE
2076  temp3 = zero
2077  END IF
2078  result( ntest ) = ( temp1+temp2 ) /
2079  $ max( unfl, temp3*ulp )
2080 *
2081  CALL clacpy( ' ', n, n, v, ldu, a, lda )
2082 *
2083 *
2084 *
2085 *
2086 * Load array V with the upper or lower triangular part
2087 * of the matrix in band form.
2088 *
2089  1190 CONTINUE
2090 *
2091  1200 CONTINUE
2092 *
2093 * End of Loop -- Check for RESULT(j) > THRESH
2094 *
2095  ntestt = ntestt + ntest
2096  CALL slafts( 'CST', n, n, jtype, ntest, result, ioldsd,
2097  $ thresh, nounit, nerrs )
2098 *
2099  1210 CONTINUE
2100  1220 CONTINUE
2101 *
2102 * Summary
2103 *
2104  CALL alasvm( 'CST', nounit, nerrs, ntestt, 0 )
2105 *
2106  9999 FORMAT( ' CDRVST2STG: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2107  $ ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
2108  9998 FORMAT( ' CDRVST2STG: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2109  $ ', KD=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
2110  $ ')' )
2111 *
2112  RETURN
2113 *
2114 * End of CDRVST2STG
2115 *
subroutine clatmr(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
CLATMR
Definition: clatmr.f:492
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine chet22(ITYPE, UPLO, N, M, KBAND, A, LDA, D, E, U, LDU, V, LDV, TAU, WORK, RWORK, RESULT)
CHET22
Definition: chet22.f:161
subroutine cheevr_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHEEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE ma...
subroutine cheevx_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
CHEEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE ma...
subroutine cheevd_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHEEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE ma...
subroutine chpevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: chpevd.f:202
subroutine cheev(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
CHEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: cheev.f:142
subroutine chbevd(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHBEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: chbevd.f:217
subroutine chpev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, RWORK, INFO)
CHPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: chpev.f:140
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
subroutine slafts(TYPE, M, N, IMAT, NTESTS, RESULT, ISEED, THRESH, IOUNIT, IE)
SLAFTS
Definition: slafts.f:101
subroutine chpevx(JOBZ, RANGE, UPLO, N, AP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)
CHPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: chpevx.f:242
subroutine chbev(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, RWORK, INFO)
CHBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: chbev.f:154
subroutine chet21(ITYPE, UPLO, N, KBAND, A, LDA, D, E, U, LDU, V, LDV, TAU, WORK, RWORK, RESULT)
CHET21
Definition: chet21.f:213
subroutine cheevr(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: cheevr.f:359
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function slarnd(IDIST, ISEED)
SLARND
Definition: slarnd.f:75
subroutine cheev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
CHEEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat...
Definition: cheev_2stage.f:191
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
real function ssxt1(IJOB, D1, N1, D2, N2, ABSTOL, ULP, UNFL)
SSXT1
Definition: ssxt1.f:108
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine cheevd(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHEEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: cheevd.f:207
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:76
subroutine chetrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
CHETRD_2STAGE
subroutine cheevx(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
CHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: cheevx.f:261
subroutine chbevd_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER...
subroutine chbevx(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)
CHBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: chbevx.f:269
subroutine chbevx_2stage(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
CHBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER...
subroutine chbev_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, RWORK, INFO)
CHBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
Definition: chbev_2stage.f:213

Here is the call graph for this function:

Here is the caller graph for this function: