LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ cdrvst2stg()

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
June 2017

Definition at line 340 of file cdrvst2stg.f.

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