LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvst ( 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 
)

CDRVST

Purpose:
      CDRVST  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 CDRVST 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,
          CDRVST 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, CDRVST
          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 CDRVST 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
November 2011

Definition at line 340 of file cdrvst.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: