LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zchkst ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
integer  NOUNIT,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( * )  AP,
double precision, dimension( * )  SD,
double precision, dimension( * )  SE,
double precision, dimension( * )  D1,
double precision, dimension( * )  D2,
double precision, dimension( * )  D3,
double precision, dimension( * )  D4,
double precision, dimension( * )  D5,
double precision, dimension( * )  WA1,
double precision, dimension( * )  WA2,
double precision, dimension( * )  WA3,
double precision, dimension( * )  WR,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( ldu, * )  V,
complex*16, dimension( * )  VP,
complex*16, dimension( * )  TAU,
complex*16, dimension( ldu, * )  Z,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

ZCHKST

Purpose:
 ZCHKST  checks the Hermitian eigenvalue problem routines.

    ZHETRD factors A as  U S U* , where * means conjugate transpose,
    S is real symmetric tridiagonal, and U is unitary.
    ZHETRD can use either just the lower or just the upper triangle
    of A; ZCHKST checks both cases.
    U is represented as a product of Householder
    transformations, whose vectors are stored in the first
    n-1 columns of V, and whose scale factors are in TAU.

    ZHPTRD does the same as ZHETRD, except that A and V are stored
    in "packed" format.

    ZUNGTR constructs the matrix U from the contents of V and TAU.

    ZUPGTR constructs the matrix U from the contents of VP and TAU.

    ZSTEQR factors S as  Z D1 Z* , where Z is the unitary
    matrix of eigenvectors and D1 is a diagonal matrix with
    the eigenvalues on the diagonal.  D2 is the matrix of
    eigenvalues computed when Z is not computed.

    DSTERF computes D3, the matrix of eigenvalues, by the
    PWK method, which does not yield eigenvectors.

    ZPTEQR factors S as  Z4 D4 Z4* , for a
    Hermitian positive definite tridiagonal matrix.
    D5 is the matrix of eigenvalues computed when Z is not
    computed.

    DSTEBZ computes selected eigenvalues.  WA1, WA2, and
    WA3 will denote eigenvalues computed to high
    absolute accuracy, with different range options.
    WR will denote eigenvalues computed to high relative
    accuracy.

    ZSTEIN computes Y, the eigenvectors of S, given the
    eigenvalues.

    ZSTEDC factors S as Z D1 Z* , where Z is the unitary
    matrix of eigenvectors and D1 is a diagonal matrix with
    the eigenvalues on the diagonal ('I' option). It may also
    update an input unitary matrix, usually the output
    from ZHETRD/ZUNGTR or ZHPTRD/ZUPGTR ('V' option). It may
    also just compute eigenvalues ('N' option).

    ZSTEMR factors S as Z D1 Z* , where Z is the unitary
    matrix of eigenvectors and D1 is a diagonal matrix with
    the eigenvalues on the diagonal ('I' option).  ZSTEMR
    uses the Relatively Robust Representation whenever possible.

 When ZCHKST 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 Hermitian eigenroutines.  For each matrix, a number
 of tests will be performed:

 (1)     | A - V S V* | / ( |A| n ulp ) ZHETRD( UPLO='U', ... )

 (2)     | I - UV* | / ( n ulp )        ZUNGTR( UPLO='U', ... )

 (3)     | A - V S V* | / ( |A| n ulp ) ZHETRD( UPLO='L', ... )

 (4)     | I - UV* | / ( n ulp )        ZUNGTR( UPLO='L', ... )

 (5-8)   Same as 1-4, but for ZHPTRD and ZUPGTR.

 (9)     | S - Z D Z* | / ( |S| n ulp ) ZSTEQR('V',...)

 (10)    | I - ZZ* | / ( n ulp )        ZSTEQR('V',...)

 (11)    | D1 - D2 | / ( |D1| ulp )        ZSTEQR('N',...)

 (12)    | D1 - D3 | / ( |D1| ulp )        DSTERF

 (13)    0 if the true eigenvalues (computed by sturm count)
         of S are within THRESH of
         those in D1.  2*THRESH if they are not.  (Tested using
         DSTECH)

 For S positive definite,

 (14)    | S - Z4 D4 Z4* | / ( |S| n ulp ) ZPTEQR('V',...)

 (15)    | I - Z4 Z4* | / ( n ulp )        ZPTEQR('V',...)

 (16)    | D4 - D5 | / ( 100 |D4| ulp )       ZPTEQR('N',...)

 When S is also diagonally dominant by the factor gamma < 1,

 (17)    max | D4(i) - WR(i) | / ( |D4(i)| omega ) ,
          i
         omega = 2 (2n-1) ULP (1 + 8 gamma**2) / (1 - gamma)**4
                                              DSTEBZ( 'A', 'E', ...)

 (18)    | WA1 - D3 | / ( |D3| ulp )          DSTEBZ( 'A', 'E', ...)

 (19)    ( max { min | WA2(i)-WA3(j) | } +
            i     j
           max { min | WA3(i)-WA2(j) | } ) / ( |D3| ulp )
            i     j
                                              DSTEBZ( 'I', 'E', ...)

 (20)    | S - Y WA1 Y* | / ( |S| n ulp )  DSTEBZ, ZSTEIN

 (21)    | I - Y Y* | / ( n ulp )          DSTEBZ, ZSTEIN

 (22)    | S - Z D Z* | / ( |S| n ulp )    ZSTEDC('I')

 (23)    | I - ZZ* | / ( n ulp )           ZSTEDC('I')

 (24)    | S - Z D Z* | / ( |S| n ulp )    ZSTEDC('V')

 (25)    | I - ZZ* | / ( n ulp )           ZSTEDC('V')

 (26)    | D1 - D2 | / ( |D1| ulp )           ZSTEDC('V') and
                                              ZSTEDC('N')

 Test 27 is disabled at the moment because ZSTEMR does not
 guarantee high relatvie accuracy.

 (27)    max | D6(i) - WR(i) | / ( |D6(i)| omega ) ,
          i
         omega = 2 (2n-1) ULP (1 + 8 gamma**2) / (1 - gamma)**4
                                              ZSTEMR('V', 'A')

 (28)    max | D6(i) - WR(i) | / ( |D6(i)| omega ) ,
          i
         omega = 2 (2n-1) ULP (1 + 8 gamma**2) / (1 - gamma)**4
                                              ZSTEMR('V', 'I')

 Tests 29 through 34 are disable at present because ZSTEMR
 does not handle partial specturm requests.

 (29)    | S - Z D Z* | / ( |S| n ulp )    ZSTEMR('V', 'I')

 (30)    | I - ZZ* | / ( n ulp )           ZSTEMR('V', 'I')

 (31)    ( max { min | WA2(i)-WA3(j) | } +
            i     j
           max { min | WA3(i)-WA2(j) | } ) / ( |D3| ulp )
            i     j
         ZSTEMR('N', 'I') vs. CSTEMR('V', 'I')

 (32)    | S - Z D Z* | / ( |S| n ulp )    ZSTEMR('V', 'V')

 (33)    | I - ZZ* | / ( n ulp )           ZSTEMR('V', 'V')

 (34)    ( max { min | WA2(i)-WA3(j) | } +
            i     j
           max { min | WA3(i)-WA2(j) | } ) / ( |D3| ulp )
            i     j
         ZSTEMR('N', 'V') vs. CSTEMR('V', 'V')

 (35)    | S - Z D Z* | / ( |S| n ulp )    ZSTEMR('V', 'A')

 (36)    | I - ZZ* | / ( n ulp )           ZSTEMR('V', 'A')

 (37)    ( max { min | WA2(i)-WA3(j) | } +
            i     j
           max { min | WA3(i)-WA2(j) | } ) / ( |D3| ulp )
            i     j
         ZSTEMR('N', 'A') vs. CSTEMR('V', 'A')

 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) Hermitian 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) Same as (8), but diagonal elements are all positive.
 (17) Same as (9), but diagonal elements are all positive.
 (18) Same as (10), but diagonal elements are all positive.
 (19) Same as (16), but multiplied by SQRT( overflow threshold )
 (20) Same as (16), but multiplied by SQRT( underflow threshold )
 (21) A diagonally dominant tridiagonal matrix with geometrically
      spaced diagonal entries 1, ..., ULP.
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of sizes of matrices to use.  If it is zero,
          ZCHKST does nothing.  It must be at least zero.
[in]NN
          NN is 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.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, ZCHKST
          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. .
[in]DOTYPE
          DOTYPE is 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.
[in,out]ISEED
          ISEED is 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 ZCHKST to continue the same random number
          sequence.
[in]THRESH
          THRESH is DOUBLE PRECISION
          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.
[in]NOUNIT
          NOUNIT is INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
[in,out]A
          A is COMPLEX*16 array of
                                  dimension ( LDA , max(NN) )
          Used to hold the matrix whose eigenvalues are to be
          computed.  On exit, A contains the last matrix actually
          used.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  It must be at
          least 1 and at least max( NN ).
[out]AP
          AP is COMPLEX*16 array of
                      dimension( max(NN)*max(NN+1)/2 )
          The matrix A stored in packed format.
[out]SD
          SD is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The diagonal of the tridiagonal matrix computed by ZHETRD.
          On exit, SD and SE contain the tridiagonal form of the
          matrix in A.
[out]SE
          SE is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The off-diagonal of the tridiagonal matrix computed by
          ZHETRD.  On exit, SD and SE contain the tridiagonal form of
          the matrix in A.
[out]D1
          D1 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The eigenvalues of A, as computed by ZSTEQR simlutaneously
          with Z.  On exit, the eigenvalues in D1 correspond with the
          matrix in A.
[out]D2
          D2 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The eigenvalues of A, as computed by ZSTEQR if Z is not
          computed.  On exit, the eigenvalues in D2 correspond with
          the matrix in A.
[out]D3
          D3 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The eigenvalues of A, as computed by DSTERF.  On exit, the
          eigenvalues in D3 correspond with the matrix in A.
[out]D4
          D4 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The eigenvalues of A, as computed by ZPTEQR(V).
          ZPTEQR factors S as  Z4 D4 Z4*
          On exit, the eigenvalues in D4 correspond with the matrix in A.
[out]D5
          D5 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The eigenvalues of A, as computed by ZPTEQR(N)
          when Z is not computed. On exit, the
          eigenvalues in D4 correspond with the matrix in A.
[out]WA1
          WA1 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          All eigenvalues of A, computed to high
          absolute accuracy, with different range options.
          as computed by DSTEBZ.
[out]WA2
          WA2 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          Selected eigenvalues of A, computed to high
          absolute accuracy, with different range options.
          as computed by DSTEBZ.
          Choose random values for IL and IU, and ask for the
          IL-th through IU-th eigenvalues.
[out]WA3
          WA3 is DOUBLE PRECISION array of
                             dimension( max(NN) )
          Selected eigenvalues of A, computed to high
          absolute accuracy, with different range options.
          as computed by DSTEBZ.
          Determine the values VL and VU of the IL-th and IU-th
          eigenvalues and ask for all eigenvalues in this range.
[out]WR
          WR is DOUBLE PRECISION array of
                             dimension( max(NN) )
          All eigenvalues of A, computed to high
          absolute accuracy, with different options.
          as computed by DSTEBZ.
[out]U
          U is COMPLEX*16 array of
                             dimension( LDU, max(NN) ).
          The unitary matrix computed by ZHETRD + ZUNGTR.
[in]LDU
          LDU is INTEGER
          The leading dimension of U, Z, and V.  It must be at least 1
          and at least max( NN ).
[out]V
          V is COMPLEX*16 array of
                             dimension( LDU, max(NN) ).
          The Housholder vectors computed by ZHETRD in reducing A to
          tridiagonal form.  The vectors computed with UPLO='U' are
          in the upper triangle, and the vectors computed with UPLO='L'
          are in the lower triangle.  (As described in ZHETRD, the
          sub- and superdiagonal are not set to 1, although the
          true Householder vector has a 1 in that position.  The
          routines that use V, such as ZUNGTR, set those entries to
          1 before using them, and then restore them later.)
[out]VP
          VP is COMPLEX*16 array of
                      dimension( max(NN)*max(NN+1)/2 )
          The matrix V stored in packed format.
[out]TAU
          TAU is COMPLEX*16 array of
                             dimension( max(NN) )
          The Householder factors computed by ZHETRD in reducing A
          to tridiagonal form.
[out]Z
          Z is COMPLEX*16 array of
                             dimension( LDU, max(NN) ).
          The unitary matrix of eigenvectors computed by ZSTEQR,
          ZPTEQR, and ZSTEIN.
[out]WORK
          WORK is COMPLEX*16 array of
                      dimension( LWORK )
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          1 + 4 * Nmax + 2 * Nmax * lg Nmax + 3 * Nmax**2
          where Nmax = max( NN(j), 2 ) and lg = log base 2.
[out]IWORK
          IWORK is INTEGER array,
          Workspace.
[out]LIWORK
          LIWORK is INTEGER
          The number of entries in IWORK.  This must be at least
                  6 + 6*Nmax + 5 * Nmax * lg Nmax 
          where Nmax = max( NN(j), 2 ) and lg = log base 2.
[out]RWORK
          RWORK is DOUBLE PRECISION array
[in]LRWORK
          LRWORK is INTEGER
          The number of entries in LRWORK (dimension( ??? )
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (26)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid
          overflow.
[out]INFO
          INFO is 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) ).
          -23: LDU < 1 or LDU < NMAX.
          -29: LWORK too small.
          If  ZLATMR, CLATMS, ZHETRD, ZUNGTR, ZSTEQR, DSTERF,
              or ZUNMC2 returns an error code, the
              absolute value of it is returned.

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

       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.
       NBLOCK          Blocksize as returned by ENVIR.
       NMAX            Largest value in NN.
       NMATS           The number of matrices generated so far.
       NERRS           The number of tests which have exceeded THRESH
                       so far.
       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 606 of file zchkst.f.

606 *
607 * -- LAPACK test routine (version 3.4.0) --
608 * -- LAPACK is a software package provided by Univ. of Tennessee, --
609 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
610 * November 2011
611 *
612 * .. Scalar Arguments ..
613  INTEGER info, lda, ldu, liwork, lrwork, lwork, nounit,
614  $ nsizes, ntypes
615  DOUBLE PRECISION thresh
616 * ..
617 * .. Array Arguments ..
618  LOGICAL dotype( * )
619  INTEGER iseed( 4 ), iwork( * ), nn( * )
620  DOUBLE PRECISION d1( * ), d2( * ), d3( * ), d4( * ), d5( * ),
621  $ result( * ), rwork( * ), sd( * ), se( * ),
622  $ wa1( * ), wa2( * ), wa3( * ), wr( * )
623  COMPLEX*16 a( lda, * ), ap( * ), tau( * ), u( ldu, * ),
624  $ v( ldu, * ), vp( * ), work( * ), z( ldu, * )
625 * ..
626 *
627 * =====================================================================
628 *
629 * .. Parameters ..
630  DOUBLE PRECISION zero, one, two, eight, ten, hun
631  parameter ( zero = 0.0d0, one = 1.0d0, two = 2.0d0,
632  $ eight = 8.0d0, ten = 10.0d0, hun = 100.0d0 )
633  COMPLEX*16 czero, cone
634  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
635  $ cone = ( 1.0d+0, 0.0d+0 ) )
636  DOUBLE PRECISION half
637  parameter ( half = one / two )
638  INTEGER maxtyp
639  parameter ( maxtyp = 21 )
640  LOGICAL crange
641  parameter ( crange = .false. )
642  LOGICAL crel
643  parameter ( crel = .false. )
644 * ..
645 * .. Local Scalars ..
646  LOGICAL badnn, tryrac
647  INTEGER i, iinfo, il, imode, inde, indrwk, itemp,
648  $ itype, iu, j, jc, jr, jsize, jtype, lgn,
649  $ liwedc, log2ui, lrwedc, lwedc, m, m2, m3,
650  $ mtypes, n, nap, nblock, nerrs, nmats, nmax,
651  $ nsplit, ntest, ntestt
652  DOUBLE PRECISION abstol, aninv, anorm, cond, ovfl, rtovfl,
653  $ rtunfl, temp1, temp2, temp3, temp4, ulp,
654  $ ulpinv, unfl, vl, vu
655 * ..
656 * .. Local Arrays ..
657  INTEGER idumma( 1 ), ioldsd( 4 ), iseed2( 4 ),
658  $ kmagn( maxtyp ), kmode( maxtyp ),
659  $ ktype( maxtyp )
660  DOUBLE PRECISION dumma( 1 )
661 * ..
662 * .. External Functions ..
663  INTEGER ilaenv
664  DOUBLE PRECISION dlamch, dlarnd, dsxt1
665  EXTERNAL ilaenv, dlamch, dlarnd, dsxt1
666 * ..
667 * .. External Subroutines ..
668  EXTERNAL dcopy, dlabad, dlasum, dstebz, dstech, dsterf,
672  $ zupgtr
673 * ..
674 * .. Intrinsic Functions ..
675  INTRINSIC abs, dble, dconjg, int, log, max, min, sqrt
676 * ..
677 * .. Data statements ..
678  DATA ktype / 1, 2, 4, 4, 4, 4, 4, 5, 5, 5, 5, 5, 8,
679  $ 8, 8, 9, 9, 9, 9, 9, 10 /
680  DATA kmagn / 1, 1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
681  $ 2, 3, 1, 1, 1, 2, 3, 1 /
682  DATA kmode / 0, 0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
683  $ 0, 0, 4, 3, 1, 4, 4, 3 /
684 * ..
685 * .. Executable Statements ..
686 *
687 * Keep ftnchek happy
688  idumma( 1 ) = 1
689 *
690 * Check for errors
691 *
692  ntestt = 0
693  info = 0
694 *
695 * Important constants
696 *
697  badnn = .false.
698  tryrac = .true.
699  nmax = 1
700  DO 10 j = 1, nsizes
701  nmax = max( nmax, nn( j ) )
702  IF( nn( j ).LT.0 )
703  $ badnn = .true.
704  10 CONTINUE
705 *
706  nblock = ilaenv( 1, 'ZHETRD', 'L', nmax, -1, -1, -1 )
707  nblock = min( nmax, max( 1, nblock ) )
708 *
709 * Check for errors
710 *
711  IF( nsizes.LT.0 ) THEN
712  info = -1
713  ELSE IF( badnn ) THEN
714  info = -2
715  ELSE IF( ntypes.LT.0 ) THEN
716  info = -3
717  ELSE IF( lda.LT.nmax ) THEN
718  info = -9
719  ELSE IF( ldu.LT.nmax ) THEN
720  info = -23
721  ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
722  info = -29
723  END IF
724 *
725  IF( info.NE.0 ) THEN
726  CALL xerbla( 'ZCHKST', -info )
727  RETURN
728  END IF
729 *
730 * Quick return if possible
731 *
732  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
733  $ RETURN
734 *
735 * More Important constants
736 *
737  unfl = dlamch( 'Safe minimum' )
738  ovfl = one / unfl
739  CALL dlabad( unfl, ovfl )
740  ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
741  ulpinv = one / ulp
742  log2ui = int( log( ulpinv ) / log( two ) )
743  rtunfl = sqrt( unfl )
744  rtovfl = sqrt( ovfl )
745 *
746 * Loop over sizes, types
747 *
748  DO 20 i = 1, 4
749  iseed2( i ) = iseed( i )
750  20 CONTINUE
751  nerrs = 0
752  nmats = 0
753 *
754  DO 310 jsize = 1, nsizes
755  n = nn( jsize )
756  IF( n.GT.0 ) THEN
757  lgn = int( log( dble( n ) ) / log( two ) )
758  IF( 2**lgn.LT.n )
759  $ lgn = lgn + 1
760  IF( 2**lgn.LT.n )
761  $ lgn = lgn + 1
762  lwedc = 1 + 4*n + 2*n*lgn + 4*n**2
763  lrwedc = 1 + 3*n + 2*n*lgn + 4*n**2
764  liwedc = 6 + 6*n + 5*n*lgn
765  ELSE
766  lwedc = 8
767  lrwedc = 7
768  liwedc = 12
769  END IF
770  nap = ( n*( n+1 ) ) / 2
771  aninv = one / dble( max( 1, n ) )
772 *
773  IF( nsizes.NE.1 ) THEN
774  mtypes = min( maxtyp, ntypes )
775  ELSE
776  mtypes = min( maxtyp+1, ntypes )
777  END IF
778 *
779  DO 300 jtype = 1, mtypes
780  IF( .NOT.dotype( jtype ) )
781  $ GO TO 300
782  nmats = nmats + 1
783  ntest = 0
784 *
785  DO 30 j = 1, 4
786  ioldsd( j ) = iseed( j )
787  30 CONTINUE
788 *
789 * Compute "A"
790 *
791 * Control parameters:
792 *
793 * KMAGN KMODE KTYPE
794 * =1 O(1) clustered 1 zero
795 * =2 large clustered 2 identity
796 * =3 small exponential (none)
797 * =4 arithmetic diagonal, (w/ eigenvalues)
798 * =5 random log Hermitian, w/ eigenvalues
799 * =6 random (none)
800 * =7 random diagonal
801 * =8 random Hermitian
802 * =9 positive definite
803 * =10 diagonally dominant tridiagonal
804 *
805  IF( mtypes.GT.maxtyp )
806  $ GO TO 100
807 *
808  itype = ktype( jtype )
809  imode = kmode( jtype )
810 *
811 * Compute norm
812 *
813  GO TO ( 40, 50, 60 )kmagn( jtype )
814 *
815  40 CONTINUE
816  anorm = one
817  GO TO 70
818 *
819  50 CONTINUE
820  anorm = ( rtovfl*ulp )*aninv
821  GO TO 70
822 *
823  60 CONTINUE
824  anorm = rtunfl*n*ulpinv
825  GO TO 70
826 *
827  70 CONTINUE
828 *
829  CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
830  iinfo = 0
831  IF( jtype.LE.15 ) THEN
832  cond = ulpinv
833  ELSE
834  cond = ulpinv*aninv / ten
835  END IF
836 *
837 * Special Matrices -- Identity & Jordan block
838 *
839 * Zero
840 *
841  IF( itype.EQ.1 ) THEN
842  iinfo = 0
843 *
844  ELSE IF( itype.EQ.2 ) THEN
845 *
846 * Identity
847 *
848  DO 80 jc = 1, n
849  a( jc, jc ) = anorm
850  80 CONTINUE
851 *
852  ELSE IF( itype.EQ.4 ) THEN
853 *
854 * Diagonal Matrix, [Eigen]values Specified
855 *
856  CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
857  $ anorm, 0, 0, 'N', a, lda, work, iinfo )
858 *
859 *
860  ELSE IF( itype.EQ.5 ) THEN
861 *
862 * Hermitian, eigenvalues specified
863 *
864  CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
865  $ anorm, n, n, 'N', a, lda, work, iinfo )
866 *
867  ELSE IF( itype.EQ.7 ) THEN
868 *
869 * Diagonal, random eigenvalues
870 *
871  CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
872  $ 'T', 'N', work( n+1 ), 1, one,
873  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
874  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
875 *
876  ELSE IF( itype.EQ.8 ) THEN
877 *
878 * Hermitian, random eigenvalues
879 *
880  CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
881  $ 'T', 'N', work( n+1 ), 1, one,
882  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
883  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
884 *
885  ELSE IF( itype.EQ.9 ) THEN
886 *
887 * Positive definite, eigenvalues specified.
888 *
889  CALL zlatms( n, n, 'S', iseed, 'P', rwork, imode, cond,
890  $ anorm, n, n, 'N', a, lda, work, iinfo )
891 *
892  ELSE IF( itype.EQ.10 ) THEN
893 *
894 * Positive definite tridiagonal, eigenvalues specified.
895 *
896  CALL zlatms( n, n, 'S', iseed, 'P', rwork, imode, cond,
897  $ anorm, 1, 1, 'N', a, lda, work, iinfo )
898  DO 90 i = 2, n
899  temp1 = abs( a( i-1, i ) )
900  temp2 = sqrt( abs( a( i-1, i-1 )*a( i, i ) ) )
901  IF( temp1.GT.half*temp2 ) THEN
902  a( i-1, i ) = a( i-1, i )*
903  $ ( half*temp2 / ( unfl+temp1 ) )
904  a( i, i-1 ) = dconjg( a( i-1, i ) )
905  END IF
906  90 CONTINUE
907 *
908  ELSE
909 *
910  iinfo = 1
911  END IF
912 *
913  IF( iinfo.NE.0 ) THEN
914  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
915  $ ioldsd
916  info = abs( iinfo )
917  RETURN
918  END IF
919 *
920  100 CONTINUE
921 *
922 * Call ZHETRD and ZUNGTR to compute S and U from
923 * upper triangle.
924 *
925  CALL zlacpy( 'U', n, n, a, lda, v, ldu )
926 *
927  ntest = 1
928  CALL zhetrd( 'U', n, v, ldu, sd, se, tau, work, lwork,
929  $ iinfo )
930 *
931  IF( iinfo.NE.0 ) THEN
932  WRITE( nounit, fmt = 9999 )'ZHETRD(U)', iinfo, n, jtype,
933  $ ioldsd
934  info = abs( iinfo )
935  IF( iinfo.LT.0 ) THEN
936  RETURN
937  ELSE
938  result( 1 ) = ulpinv
939  GO TO 280
940  END IF
941  END IF
942 *
943  CALL zlacpy( 'U', n, n, v, ldu, u, ldu )
944 *
945  ntest = 2
946  CALL zungtr( 'U', n, u, ldu, tau, work, lwork, iinfo )
947  IF( iinfo.NE.0 ) THEN
948  WRITE( nounit, fmt = 9999 )'ZUNGTR(U)', iinfo, n, jtype,
949  $ ioldsd
950  info = abs( iinfo )
951  IF( iinfo.LT.0 ) THEN
952  RETURN
953  ELSE
954  result( 2 ) = ulpinv
955  GO TO 280
956  END IF
957  END IF
958 *
959 * Do tests 1 and 2
960 *
961  CALL zhet21( 2, 'Upper', n, 1, a, lda, sd, se, u, ldu, v,
962  $ ldu, tau, work, rwork, result( 1 ) )
963  CALL zhet21( 3, 'Upper', n, 1, a, lda, sd, se, u, ldu, v,
964  $ ldu, tau, work, rwork, result( 2 ) )
965 *
966 * Call ZHETRD and ZUNGTR to compute S and U from
967 * lower triangle, do tests.
968 *
969  CALL zlacpy( 'L', n, n, a, lda, v, ldu )
970 *
971  ntest = 3
972  CALL zhetrd( 'L', n, v, ldu, sd, se, tau, work, lwork,
973  $ iinfo )
974 *
975  IF( iinfo.NE.0 ) THEN
976  WRITE( nounit, fmt = 9999 )'ZHETRD(L)', iinfo, n, jtype,
977  $ ioldsd
978  info = abs( iinfo )
979  IF( iinfo.LT.0 ) THEN
980  RETURN
981  ELSE
982  result( 3 ) = ulpinv
983  GO TO 280
984  END IF
985  END IF
986 *
987  CALL zlacpy( 'L', n, n, v, ldu, u, ldu )
988 *
989  ntest = 4
990  CALL zungtr( 'L', n, u, ldu, tau, work, lwork, iinfo )
991  IF( iinfo.NE.0 ) THEN
992  WRITE( nounit, fmt = 9999 )'ZUNGTR(L)', iinfo, n, jtype,
993  $ ioldsd
994  info = abs( iinfo )
995  IF( iinfo.LT.0 ) THEN
996  RETURN
997  ELSE
998  result( 4 ) = ulpinv
999  GO TO 280
1000  END IF
1001  END IF
1002 *
1003  CALL zhet21( 2, 'Lower', n, 1, a, lda, sd, se, u, ldu, v,
1004  $ ldu, tau, work, rwork, result( 3 ) )
1005  CALL zhet21( 3, 'Lower', n, 1, a, lda, sd, se, u, ldu, v,
1006  $ ldu, tau, work, rwork, result( 4 ) )
1007 *
1008 * Store the upper triangle of A in AP
1009 *
1010  i = 0
1011  DO 120 jc = 1, n
1012  DO 110 jr = 1, jc
1013  i = i + 1
1014  ap( i ) = a( jr, jc )
1015  110 CONTINUE
1016  120 CONTINUE
1017 *
1018 * Call ZHPTRD and ZUPGTR to compute S and U from AP
1019 *
1020  CALL zcopy( nap, ap, 1, vp, 1 )
1021 *
1022  ntest = 5
1023  CALL zhptrd( 'U', n, vp, sd, se, tau, iinfo )
1024 *
1025  IF( iinfo.NE.0 ) THEN
1026  WRITE( nounit, fmt = 9999 )'ZHPTRD(U)', iinfo, n, jtype,
1027  $ ioldsd
1028  info = abs( iinfo )
1029  IF( iinfo.LT.0 ) THEN
1030  RETURN
1031  ELSE
1032  result( 5 ) = ulpinv
1033  GO TO 280
1034  END IF
1035  END IF
1036 *
1037  ntest = 6
1038  CALL zupgtr( 'U', n, vp, tau, u, ldu, work, iinfo )
1039  IF( iinfo.NE.0 ) THEN
1040  WRITE( nounit, fmt = 9999 )'ZUPGTR(U)', iinfo, n, jtype,
1041  $ ioldsd
1042  info = abs( iinfo )
1043  IF( iinfo.LT.0 ) THEN
1044  RETURN
1045  ELSE
1046  result( 6 ) = ulpinv
1047  GO TO 280
1048  END IF
1049  END IF
1050 *
1051 * Do tests 5 and 6
1052 *
1053  CALL zhpt21( 2, 'Upper', n, 1, ap, sd, se, u, ldu, vp, tau,
1054  $ work, rwork, result( 5 ) )
1055  CALL zhpt21( 3, 'Upper', n, 1, ap, sd, se, u, ldu, vp, tau,
1056  $ work, rwork, result( 6 ) )
1057 *
1058 * Store the lower triangle of A in AP
1059 *
1060  i = 0
1061  DO 140 jc = 1, n
1062  DO 130 jr = jc, n
1063  i = i + 1
1064  ap( i ) = a( jr, jc )
1065  130 CONTINUE
1066  140 CONTINUE
1067 *
1068 * Call ZHPTRD and ZUPGTR to compute S and U from AP
1069 *
1070  CALL zcopy( nap, ap, 1, vp, 1 )
1071 *
1072  ntest = 7
1073  CALL zhptrd( 'L', n, vp, sd, se, tau, iinfo )
1074 *
1075  IF( iinfo.NE.0 ) THEN
1076  WRITE( nounit, fmt = 9999 )'ZHPTRD(L)', iinfo, n, jtype,
1077  $ ioldsd
1078  info = abs( iinfo )
1079  IF( iinfo.LT.0 ) THEN
1080  RETURN
1081  ELSE
1082  result( 7 ) = ulpinv
1083  GO TO 280
1084  END IF
1085  END IF
1086 *
1087  ntest = 8
1088  CALL zupgtr( 'L', n, vp, tau, u, ldu, work, iinfo )
1089  IF( iinfo.NE.0 ) THEN
1090  WRITE( nounit, fmt = 9999 )'ZUPGTR(L)', iinfo, n, jtype,
1091  $ ioldsd
1092  info = abs( iinfo )
1093  IF( iinfo.LT.0 ) THEN
1094  RETURN
1095  ELSE
1096  result( 8 ) = ulpinv
1097  GO TO 280
1098  END IF
1099  END IF
1100 *
1101  CALL zhpt21( 2, 'Lower', n, 1, ap, sd, se, u, ldu, vp, tau,
1102  $ work, rwork, result( 7 ) )
1103  CALL zhpt21( 3, 'Lower', n, 1, ap, sd, se, u, ldu, vp, tau,
1104  $ work, rwork, result( 8 ) )
1105 *
1106 * Call ZSTEQR to compute D1, D2, and Z, do tests.
1107 *
1108 * Compute D1 and Z
1109 *
1110  CALL dcopy( n, sd, 1, d1, 1 )
1111  IF( n.GT.0 )
1112  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1113  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1114 *
1115  ntest = 9
1116  CALL zsteqr( 'V', n, d1, rwork, z, ldu, rwork( n+1 ),
1117  $ iinfo )
1118  IF( iinfo.NE.0 ) THEN
1119  WRITE( nounit, fmt = 9999 )'ZSTEQR(V)', iinfo, n, jtype,
1120  $ ioldsd
1121  info = abs( iinfo )
1122  IF( iinfo.LT.0 ) THEN
1123  RETURN
1124  ELSE
1125  result( 9 ) = ulpinv
1126  GO TO 280
1127  END IF
1128  END IF
1129 *
1130 * Compute D2
1131 *
1132  CALL dcopy( n, sd, 1, d2, 1 )
1133  IF( n.GT.0 )
1134  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1135 *
1136  ntest = 11
1137  CALL zsteqr( 'N', n, d2, rwork, work, ldu, rwork( n+1 ),
1138  $ iinfo )
1139  IF( iinfo.NE.0 ) THEN
1140  WRITE( nounit, fmt = 9999 )'ZSTEQR(N)', iinfo, n, jtype,
1141  $ ioldsd
1142  info = abs( iinfo )
1143  IF( iinfo.LT.0 ) THEN
1144  RETURN
1145  ELSE
1146  result( 11 ) = ulpinv
1147  GO TO 280
1148  END IF
1149  END IF
1150 *
1151 * Compute D3 (using PWK method)
1152 *
1153  CALL dcopy( n, sd, 1, d3, 1 )
1154  IF( n.GT.0 )
1155  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1156 *
1157  ntest = 12
1158  CALL dsterf( n, d3, rwork, iinfo )
1159  IF( iinfo.NE.0 ) THEN
1160  WRITE( nounit, fmt = 9999 )'DSTERF', iinfo, n, jtype,
1161  $ ioldsd
1162  info = abs( iinfo )
1163  IF( iinfo.LT.0 ) THEN
1164  RETURN
1165  ELSE
1166  result( 12 ) = ulpinv
1167  GO TO 280
1168  END IF
1169  END IF
1170 *
1171 * Do Tests 9 and 10
1172 *
1173  CALL zstt21( n, 0, sd, se, d1, dumma, z, ldu, work, rwork,
1174  $ result( 9 ) )
1175 *
1176 * Do Tests 11 and 12
1177 *
1178  temp1 = zero
1179  temp2 = zero
1180  temp3 = zero
1181  temp4 = zero
1182 *
1183  DO 150 j = 1, n
1184  temp1 = max( temp1, abs( d1( j ) ), abs( d2( j ) ) )
1185  temp2 = max( temp2, abs( d1( j )-d2( j ) ) )
1186  temp3 = max( temp3, abs( d1( j ) ), abs( d3( j ) ) )
1187  temp4 = max( temp4, abs( d1( j )-d3( j ) ) )
1188  150 CONTINUE
1189 *
1190  result( 11 ) = temp2 / max( unfl, ulp*max( temp1, temp2 ) )
1191  result( 12 ) = temp4 / max( unfl, ulp*max( temp3, temp4 ) )
1192 *
1193 * Do Test 13 -- Sturm Sequence Test of Eigenvalues
1194 * Go up by factors of two until it succeeds
1195 *
1196  ntest = 13
1197  temp1 = thresh*( half-ulp )
1198 *
1199  DO 160 j = 0, log2ui
1200  CALL dstech( n, sd, se, d1, temp1, rwork, iinfo )
1201  IF( iinfo.EQ.0 )
1202  $ GO TO 170
1203  temp1 = temp1*two
1204  160 CONTINUE
1205 *
1206  170 CONTINUE
1207  result( 13 ) = temp1
1208 *
1209 * For positive definite matrices ( JTYPE.GT.15 ) call ZPTEQR
1210 * and do tests 14, 15, and 16 .
1211 *
1212  IF( jtype.GT.15 ) THEN
1213 *
1214 * Compute D4 and Z4
1215 *
1216  CALL dcopy( n, sd, 1, d4, 1 )
1217  IF( n.GT.0 )
1218  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1219  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1220 *
1221  ntest = 14
1222  CALL zpteqr( 'V', n, d4, rwork, z, ldu, rwork( n+1 ),
1223  $ iinfo )
1224  IF( iinfo.NE.0 ) THEN
1225  WRITE( nounit, fmt = 9999 )'ZPTEQR(V)', iinfo, n,
1226  $ jtype, ioldsd
1227  info = abs( iinfo )
1228  IF( iinfo.LT.0 ) THEN
1229  RETURN
1230  ELSE
1231  result( 14 ) = ulpinv
1232  GO TO 280
1233  END IF
1234  END IF
1235 *
1236 * Do Tests 14 and 15
1237 *
1238  CALL zstt21( n, 0, sd, se, d4, dumma, z, ldu, work,
1239  $ rwork, result( 14 ) )
1240 *
1241 * Compute D5
1242 *
1243  CALL dcopy( n, sd, 1, d5, 1 )
1244  IF( n.GT.0 )
1245  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1246 *
1247  ntest = 16
1248  CALL zpteqr( 'N', n, d5, rwork, z, ldu, rwork( n+1 ),
1249  $ iinfo )
1250  IF( iinfo.NE.0 ) THEN
1251  WRITE( nounit, fmt = 9999 )'ZPTEQR(N)', iinfo, n,
1252  $ jtype, ioldsd
1253  info = abs( iinfo )
1254  IF( iinfo.LT.0 ) THEN
1255  RETURN
1256  ELSE
1257  result( 16 ) = ulpinv
1258  GO TO 280
1259  END IF
1260  END IF
1261 *
1262 * Do Test 16
1263 *
1264  temp1 = zero
1265  temp2 = zero
1266  DO 180 j = 1, n
1267  temp1 = max( temp1, abs( d4( j ) ), abs( d5( j ) ) )
1268  temp2 = max( temp2, abs( d4( j )-d5( j ) ) )
1269  180 CONTINUE
1270 *
1271  result( 16 ) = temp2 / max( unfl,
1272  $ hun*ulp*max( temp1, temp2 ) )
1273  ELSE
1274  result( 14 ) = zero
1275  result( 15 ) = zero
1276  result( 16 ) = zero
1277  END IF
1278 *
1279 * Call DSTEBZ with different options and do tests 17-18.
1280 *
1281 * If S is positive definite and diagonally dominant,
1282 * ask for all eigenvalues with high relative accuracy.
1283 *
1284  vl = zero
1285  vu = zero
1286  il = 0
1287  iu = 0
1288  IF( jtype.EQ.21 ) THEN
1289  ntest = 17
1290  abstol = unfl + unfl
1291  CALL dstebz( 'A', 'E', n, vl, vu, il, iu, abstol, sd, se,
1292  $ m, nsplit, wr, iwork( 1 ), iwork( n+1 ),
1293  $ rwork, iwork( 2*n+1 ), iinfo )
1294  IF( iinfo.NE.0 ) THEN
1295  WRITE( nounit, fmt = 9999 )'DSTEBZ(A,rel)', iinfo, n,
1296  $ jtype, ioldsd
1297  info = abs( iinfo )
1298  IF( iinfo.LT.0 ) THEN
1299  RETURN
1300  ELSE
1301  result( 17 ) = ulpinv
1302  GO TO 280
1303  END IF
1304  END IF
1305 *
1306 * Do test 17
1307 *
1308  temp2 = two*( two*n-one )*ulp*( one+eight*half**2 ) /
1309  $ ( one-half )**4
1310 *
1311  temp1 = zero
1312  DO 190 j = 1, n
1313  temp1 = max( temp1, abs( d4( j )-wr( n-j+1 ) ) /
1314  $ ( abstol+abs( d4( j ) ) ) )
1315  190 CONTINUE
1316 *
1317  result( 17 ) = temp1 / temp2
1318  ELSE
1319  result( 17 ) = zero
1320  END IF
1321 *
1322 * Now ask for all eigenvalues with high absolute accuracy.
1323 *
1324  ntest = 18
1325  abstol = unfl + unfl
1326  CALL dstebz( 'A', 'E', n, vl, vu, il, iu, abstol, sd, se, m,
1327  $ nsplit, wa1, iwork( 1 ), iwork( n+1 ), rwork,
1328  $ iwork( 2*n+1 ), iinfo )
1329  IF( iinfo.NE.0 ) THEN
1330  WRITE( nounit, fmt = 9999 )'DSTEBZ(A)', iinfo, n, jtype,
1331  $ ioldsd
1332  info = abs( iinfo )
1333  IF( iinfo.LT.0 ) THEN
1334  RETURN
1335  ELSE
1336  result( 18 ) = ulpinv
1337  GO TO 280
1338  END IF
1339  END IF
1340 *
1341 * Do test 18
1342 *
1343  temp1 = zero
1344  temp2 = zero
1345  DO 200 j = 1, n
1346  temp1 = max( temp1, abs( d3( j ) ), abs( wa1( j ) ) )
1347  temp2 = max( temp2, abs( d3( j )-wa1( j ) ) )
1348  200 CONTINUE
1349 *
1350  result( 18 ) = temp2 / max( unfl, ulp*max( temp1, temp2 ) )
1351 *
1352 * Choose random values for IL and IU, and ask for the
1353 * IL-th through IU-th eigenvalues.
1354 *
1355  ntest = 19
1356  IF( n.LE.1 ) THEN
1357  il = 1
1358  iu = n
1359  ELSE
1360  il = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1361  iu = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1362  IF( iu.LT.il ) THEN
1363  itemp = iu
1364  iu = il
1365  il = itemp
1366  END IF
1367  END IF
1368 *
1369  CALL dstebz( 'I', 'E', n, vl, vu, il, iu, abstol, sd, se,
1370  $ m2, nsplit, wa2, iwork( 1 ), iwork( n+1 ),
1371  $ rwork, iwork( 2*n+1 ), iinfo )
1372  IF( iinfo.NE.0 ) THEN
1373  WRITE( nounit, fmt = 9999 )'DSTEBZ(I)', iinfo, n, jtype,
1374  $ ioldsd
1375  info = abs( iinfo )
1376  IF( iinfo.LT.0 ) THEN
1377  RETURN
1378  ELSE
1379  result( 19 ) = ulpinv
1380  GO TO 280
1381  END IF
1382  END IF
1383 *
1384 * Determine the values VL and VU of the IL-th and IU-th
1385 * eigenvalues and ask for all eigenvalues in this range.
1386 *
1387  IF( n.GT.0 ) THEN
1388  IF( il.NE.1 ) THEN
1389  vl = wa1( il ) - max( half*( wa1( il )-wa1( il-1 ) ),
1390  $ ulp*anorm, two*rtunfl )
1391  ELSE
1392  vl = wa1( 1 ) - max( half*( wa1( n )-wa1( 1 ) ),
1393  $ ulp*anorm, two*rtunfl )
1394  END IF
1395  IF( iu.NE.n ) THEN
1396  vu = wa1( iu ) + max( half*( wa1( iu+1 )-wa1( iu ) ),
1397  $ ulp*anorm, two*rtunfl )
1398  ELSE
1399  vu = wa1( n ) + max( half*( wa1( n )-wa1( 1 ) ),
1400  $ ulp*anorm, two*rtunfl )
1401  END IF
1402  ELSE
1403  vl = zero
1404  vu = one
1405  END IF
1406 *
1407  CALL dstebz( 'V', 'E', n, vl, vu, il, iu, abstol, sd, se,
1408  $ m3, nsplit, wa3, iwork( 1 ), iwork( n+1 ),
1409  $ rwork, iwork( 2*n+1 ), iinfo )
1410  IF( iinfo.NE.0 ) THEN
1411  WRITE( nounit, fmt = 9999 )'DSTEBZ(V)', iinfo, n, jtype,
1412  $ ioldsd
1413  info = abs( iinfo )
1414  IF( iinfo.LT.0 ) THEN
1415  RETURN
1416  ELSE
1417  result( 19 ) = ulpinv
1418  GO TO 280
1419  END IF
1420  END IF
1421 *
1422  IF( m3.EQ.0 .AND. n.NE.0 ) THEN
1423  result( 19 ) = ulpinv
1424  GO TO 280
1425  END IF
1426 *
1427 * Do test 19
1428 *
1429  temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1430  temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1431  IF( n.GT.0 ) THEN
1432  temp3 = max( abs( wa1( n ) ), abs( wa1( 1 ) ) )
1433  ELSE
1434  temp3 = zero
1435  END IF
1436 *
1437  result( 19 ) = ( temp1+temp2 ) / max( unfl, temp3*ulp )
1438 *
1439 * Call ZSTEIN to compute eigenvectors corresponding to
1440 * eigenvalues in WA1. (First call DSTEBZ again, to make sure
1441 * it returns these eigenvalues in the correct order.)
1442 *
1443  ntest = 21
1444  CALL dstebz( 'A', 'B', n, vl, vu, il, iu, abstol, sd, se, m,
1445  $ nsplit, wa1, iwork( 1 ), iwork( n+1 ), rwork,
1446  $ iwork( 2*n+1 ), iinfo )
1447  IF( iinfo.NE.0 ) THEN
1448  WRITE( nounit, fmt = 9999 )'DSTEBZ(A,B)', iinfo, n,
1449  $ jtype, ioldsd
1450  info = abs( iinfo )
1451  IF( iinfo.LT.0 ) THEN
1452  RETURN
1453  ELSE
1454  result( 20 ) = ulpinv
1455  result( 21 ) = ulpinv
1456  GO TO 280
1457  END IF
1458  END IF
1459 *
1460  CALL zstein( n, sd, se, m, wa1, iwork( 1 ), iwork( n+1 ), z,
1461  $ ldu, rwork, iwork( 2*n+1 ), iwork( 3*n+1 ),
1462  $ iinfo )
1463  IF( iinfo.NE.0 ) THEN
1464  WRITE( nounit, fmt = 9999 )'ZSTEIN', iinfo, n, jtype,
1465  $ ioldsd
1466  info = abs( iinfo )
1467  IF( iinfo.LT.0 ) THEN
1468  RETURN
1469  ELSE
1470  result( 20 ) = ulpinv
1471  result( 21 ) = ulpinv
1472  GO TO 280
1473  END IF
1474  END IF
1475 *
1476 * Do tests 20 and 21
1477 *
1478  CALL zstt21( n, 0, sd, se, wa1, dumma, z, ldu, work, rwork,
1479  $ result( 20 ) )
1480 *
1481 * Call ZSTEDC(I) to compute D1 and Z, do tests.
1482 *
1483 * Compute D1 and Z
1484 *
1485  inde = 1
1486  indrwk = inde + n
1487  CALL dcopy( n, sd, 1, d1, 1 )
1488  IF( n.GT.0 )
1489  $ CALL dcopy( n-1, se, 1, rwork( inde ), 1 )
1490  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1491 *
1492  ntest = 22
1493  CALL zstedc( 'I', n, d1, rwork( inde ), z, ldu, work, lwedc,
1494  $ rwork( indrwk ), lrwedc, iwork, liwedc, iinfo )
1495  IF( iinfo.NE.0 ) THEN
1496  WRITE( nounit, fmt = 9999 )'ZSTEDC(I)', iinfo, n, jtype,
1497  $ ioldsd
1498  info = abs( iinfo )
1499  IF( iinfo.LT.0 ) THEN
1500  RETURN
1501  ELSE
1502  result( 22 ) = ulpinv
1503  GO TO 280
1504  END IF
1505  END IF
1506 *
1507 * Do Tests 22 and 23
1508 *
1509  CALL zstt21( n, 0, sd, se, d1, dumma, z, ldu, work, rwork,
1510  $ result( 22 ) )
1511 *
1512 * Call ZSTEDC(V) to compute D1 and Z, do tests.
1513 *
1514 * Compute D1 and Z
1515 *
1516  CALL dcopy( n, sd, 1, d1, 1 )
1517  IF( n.GT.0 )
1518  $ CALL dcopy( n-1, se, 1, rwork( inde ), 1 )
1519  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1520 *
1521  ntest = 24
1522  CALL zstedc( 'V', n, d1, rwork( inde ), z, ldu, work, lwedc,
1523  $ rwork( indrwk ), lrwedc, iwork, liwedc, iinfo )
1524  IF( iinfo.NE.0 ) THEN
1525  WRITE( nounit, fmt = 9999 )'ZSTEDC(V)', iinfo, n, jtype,
1526  $ ioldsd
1527  info = abs( iinfo )
1528  IF( iinfo.LT.0 ) THEN
1529  RETURN
1530  ELSE
1531  result( 24 ) = ulpinv
1532  GO TO 280
1533  END IF
1534  END IF
1535 *
1536 * Do Tests 24 and 25
1537 *
1538  CALL zstt21( n, 0, sd, se, d1, dumma, z, ldu, work, rwork,
1539  $ result( 24 ) )
1540 *
1541 * Call ZSTEDC(N) to compute D2, do tests.
1542 *
1543 * Compute D2
1544 *
1545  CALL dcopy( n, sd, 1, d2, 1 )
1546  IF( n.GT.0 )
1547  $ CALL dcopy( n-1, se, 1, rwork( inde ), 1 )
1548  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1549 *
1550  ntest = 26
1551  CALL zstedc( 'N', n, d2, rwork( inde ), z, ldu, work, lwedc,
1552  $ rwork( indrwk ), lrwedc, iwork, liwedc, iinfo )
1553  IF( iinfo.NE.0 ) THEN
1554  WRITE( nounit, fmt = 9999 )'ZSTEDC(N)', iinfo, n, jtype,
1555  $ ioldsd
1556  info = abs( iinfo )
1557  IF( iinfo.LT.0 ) THEN
1558  RETURN
1559  ELSE
1560  result( 26 ) = ulpinv
1561  GO TO 280
1562  END IF
1563  END IF
1564 *
1565 * Do Test 26
1566 *
1567  temp1 = zero
1568  temp2 = zero
1569 *
1570  DO 210 j = 1, n
1571  temp1 = max( temp1, abs( d1( j ) ), abs( d2( j ) ) )
1572  temp2 = max( temp2, abs( d1( j )-d2( j ) ) )
1573  210 CONTINUE
1574 *
1575  result( 26 ) = temp2 / max( unfl, ulp*max( temp1, temp2 ) )
1576 *
1577 * Only test ZSTEMR if IEEE compliant
1578 *
1579  IF( ilaenv( 10, 'ZSTEMR', 'VA', 1, 0, 0, 0 ).EQ.1 .AND.
1580  $ ilaenv( 11, 'ZSTEMR', 'VA', 1, 0, 0, 0 ).EQ.1 ) THEN
1581 *
1582 * Call ZSTEMR, do test 27 (relative eigenvalue accuracy)
1583 *
1584 * If S is positive definite and diagonally dominant,
1585 * ask for all eigenvalues with high relative accuracy.
1586 *
1587  vl = zero
1588  vu = zero
1589  il = 0
1590  iu = 0
1591  IF( jtype.EQ.21 .AND. crel ) THEN
1592  ntest = 27
1593  abstol = unfl + unfl
1594  CALL zstemr( 'V', 'A', n, sd, se, vl, vu, il, iu,
1595  $ m, wr, z, ldu, n, iwork( 1 ), tryrac,
1596  $ rwork, lrwork, iwork( 2*n+1 ), lwork-2*n,
1597  $ iinfo )
1598  IF( iinfo.NE.0 ) THEN
1599  WRITE( nounit, fmt = 9999 )'ZSTEMR(V,A,rel)',
1600  $ iinfo, n, jtype, ioldsd
1601  info = abs( iinfo )
1602  IF( iinfo.LT.0 ) THEN
1603  RETURN
1604  ELSE
1605  result( 27 ) = ulpinv
1606  GO TO 270
1607  END IF
1608  END IF
1609 *
1610 * Do test 27
1611 *
1612  temp2 = two*( two*n-one )*ulp*( one+eight*half**2 ) /
1613  $ ( one-half )**4
1614 *
1615  temp1 = zero
1616  DO 220 j = 1, n
1617  temp1 = max( temp1, abs( d4( j )-wr( n-j+1 ) ) /
1618  $ ( abstol+abs( d4( j ) ) ) )
1619  220 CONTINUE
1620 *
1621  result( 27 ) = temp1 / temp2
1622 *
1623  il = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1624  iu = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1625  IF( iu.LT.il ) THEN
1626  itemp = iu
1627  iu = il
1628  il = itemp
1629  END IF
1630 *
1631  IF( crange ) THEN
1632  ntest = 28
1633  abstol = unfl + unfl
1634  CALL zstemr( 'V', 'I', n, sd, se, vl, vu, il, iu,
1635  $ m, wr, z, ldu, n, iwork( 1 ), tryrac,
1636  $ rwork, lrwork, iwork( 2*n+1 ),
1637  $ lwork-2*n, iinfo )
1638 *
1639  IF( iinfo.NE.0 ) THEN
1640  WRITE( nounit, fmt = 9999 )'ZSTEMR(V,I,rel)',
1641  $ iinfo, n, jtype, ioldsd
1642  info = abs( iinfo )
1643  IF( iinfo.LT.0 ) THEN
1644  RETURN
1645  ELSE
1646  result( 28 ) = ulpinv
1647  GO TO 270
1648  END IF
1649  END IF
1650 *
1651 *
1652 * Do test 28
1653 *
1654  temp2 = two*( two*n-one )*ulp*
1655  $ ( one+eight*half**2 ) / ( one-half )**4
1656 *
1657  temp1 = zero
1658  DO 230 j = il, iu
1659  temp1 = max( temp1, abs( wr( j-il+1 )-d4( n-j+
1660  $ 1 ) ) / ( abstol+abs( wr( j-il+1 ) ) ) )
1661  230 CONTINUE
1662 *
1663  result( 28 ) = temp1 / temp2
1664  ELSE
1665  result( 28 ) = zero
1666  END IF
1667  ELSE
1668  result( 27 ) = zero
1669  result( 28 ) = zero
1670  END IF
1671 *
1672 * Call ZSTEMR(V,I) to compute D1 and Z, do tests.
1673 *
1674 * Compute D1 and Z
1675 *
1676  CALL dcopy( n, sd, 1, d5, 1 )
1677  IF( n.GT.0 )
1678  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1679  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1680 *
1681  IF( crange ) THEN
1682  ntest = 29
1683  il = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1684  iu = 1 + ( n-1 )*int( dlarnd( 1, iseed2 ) )
1685  IF( iu.LT.il ) THEN
1686  itemp = iu
1687  iu = il
1688  il = itemp
1689  END IF
1690  CALL zstemr( 'V', 'I', n, d5, rwork, vl, vu, il, iu,
1691  $ m, d1, z, ldu, n, iwork( 1 ), tryrac,
1692  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1693  $ liwork-2*n, iinfo )
1694  IF( iinfo.NE.0 ) THEN
1695  WRITE( nounit, fmt = 9999 )'ZSTEMR(V,I)', iinfo,
1696  $ n, jtype, ioldsd
1697  info = abs( iinfo )
1698  IF( iinfo.LT.0 ) THEN
1699  RETURN
1700  ELSE
1701  result( 29 ) = ulpinv
1702  GO TO 280
1703  END IF
1704  END IF
1705 *
1706 * Do Tests 29 and 30
1707 *
1708 *
1709 * Call ZSTEMR to compute D2, do tests.
1710 *
1711 * Compute D2
1712 *
1713  CALL dcopy( n, sd, 1, d5, 1 )
1714  IF( n.GT.0 )
1715  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1716 *
1717  ntest = 31
1718  CALL zstemr( 'N', 'I', n, d5, rwork, vl, vu, il, iu,
1719  $ m, d2, z, ldu, n, iwork( 1 ), tryrac,
1720  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1721  $ liwork-2*n, iinfo )
1722  IF( iinfo.NE.0 ) THEN
1723  WRITE( nounit, fmt = 9999 )'ZSTEMR(N,I)', iinfo,
1724  $ n, jtype, ioldsd
1725  info = abs( iinfo )
1726  IF( iinfo.LT.0 ) THEN
1727  RETURN
1728  ELSE
1729  result( 31 ) = ulpinv
1730  GO TO 280
1731  END IF
1732  END IF
1733 *
1734 * Do Test 31
1735 *
1736  temp1 = zero
1737  temp2 = zero
1738 *
1739  DO 240 j = 1, iu - il + 1
1740  temp1 = max( temp1, abs( d1( j ) ),
1741  $ abs( d2( j ) ) )
1742  temp2 = max( temp2, abs( d1( j )-d2( j ) ) )
1743  240 CONTINUE
1744 *
1745  result( 31 ) = temp2 / max( unfl,
1746  $ ulp*max( temp1, temp2 ) )
1747 *
1748 *
1749 * Call ZSTEMR(V,V) to compute D1 and Z, do tests.
1750 *
1751 * Compute D1 and Z
1752 *
1753  CALL dcopy( n, sd, 1, d5, 1 )
1754  IF( n.GT.0 )
1755  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1756  CALL zlaset( 'Full', n, n, czero, cone, z, ldu )
1757 *
1758  ntest = 32
1759 *
1760  IF( n.GT.0 ) THEN
1761  IF( il.NE.1 ) THEN
1762  vl = d2( il ) - max( half*
1763  $ ( d2( il )-d2( il-1 ) ), ulp*anorm,
1764  $ two*rtunfl )
1765  ELSE
1766  vl = d2( 1 ) - max( half*( d2( n )-d2( 1 ) ),
1767  $ ulp*anorm, two*rtunfl )
1768  END IF
1769  IF( iu.NE.n ) THEN
1770  vu = d2( iu ) + max( half*
1771  $ ( d2( iu+1 )-d2( iu ) ), ulp*anorm,
1772  $ two*rtunfl )
1773  ELSE
1774  vu = d2( n ) + max( half*( d2( n )-d2( 1 ) ),
1775  $ ulp*anorm, two*rtunfl )
1776  END IF
1777  ELSE
1778  vl = zero
1779  vu = one
1780  END IF
1781 *
1782  CALL zstemr( 'V', 'V', n, d5, rwork, vl, vu, il, iu,
1783  $ m, d1, z, ldu, m, iwork( 1 ), tryrac,
1784  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1785  $ liwork-2*n, iinfo )
1786  IF( iinfo.NE.0 ) THEN
1787  WRITE( nounit, fmt = 9999 )'ZSTEMR(V,V)', iinfo,
1788  $ n, jtype, ioldsd
1789  info = abs( iinfo )
1790  IF( iinfo.LT.0 ) THEN
1791  RETURN
1792  ELSE
1793  result( 32 ) = ulpinv
1794  GO TO 280
1795  END IF
1796  END IF
1797 *
1798 * Do Tests 32 and 33
1799 *
1800  CALL zstt22( n, m, 0, sd, se, d1, dumma, z, ldu, work,
1801  $ m, rwork, result( 32 ) )
1802 *
1803 * Call ZSTEMR to compute D2, do tests.
1804 *
1805 * Compute D2
1806 *
1807  CALL dcopy( n, sd, 1, d5, 1 )
1808  IF( n.GT.0 )
1809  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1810 *
1811  ntest = 34
1812  CALL zstemr( 'N', 'V', n, d5, rwork, vl, vu, il, iu,
1813  $ m, d2, z, ldu, n, iwork( 1 ), tryrac,
1814  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1815  $ liwork-2*n, iinfo )
1816  IF( iinfo.NE.0 ) THEN
1817  WRITE( nounit, fmt = 9999 )'ZSTEMR(N,V)', iinfo,
1818  $ n, jtype, ioldsd
1819  info = abs( iinfo )
1820  IF( iinfo.LT.0 ) THEN
1821  RETURN
1822  ELSE
1823  result( 34 ) = ulpinv
1824  GO TO 280
1825  END IF
1826  END IF
1827 *
1828 * Do Test 34
1829 *
1830  temp1 = zero
1831  temp2 = zero
1832 *
1833  DO 250 j = 1, iu - il + 1
1834  temp1 = max( temp1, abs( d1( j ) ),
1835  $ abs( d2( j ) ) )
1836  temp2 = max( temp2, abs( d1( j )-d2( j ) ) )
1837  250 CONTINUE
1838 *
1839  result( 34 ) = temp2 / max( unfl,
1840  $ ulp*max( temp1, temp2 ) )
1841  ELSE
1842  result( 29 ) = zero
1843  result( 30 ) = zero
1844  result( 31 ) = zero
1845  result( 32 ) = zero
1846  result( 33 ) = zero
1847  result( 34 ) = zero
1848  END IF
1849 *
1850 *
1851 * Call ZSTEMR(V,A) to compute D1 and Z, do tests.
1852 *
1853 * Compute D1 and Z
1854 *
1855  CALL dcopy( n, sd, 1, d5, 1 )
1856  IF( n.GT.0 )
1857  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1858 *
1859  ntest = 35
1860 *
1861  CALL zstemr( 'V', 'A', n, d5, rwork, vl, vu, il, iu,
1862  $ m, d1, z, ldu, n, iwork( 1 ), tryrac,
1863  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1864  $ liwork-2*n, iinfo )
1865  IF( iinfo.NE.0 ) THEN
1866  WRITE( nounit, fmt = 9999 )'ZSTEMR(V,A)', iinfo, n,
1867  $ jtype, ioldsd
1868  info = abs( iinfo )
1869  IF( iinfo.LT.0 ) THEN
1870  RETURN
1871  ELSE
1872  result( 35 ) = ulpinv
1873  GO TO 280
1874  END IF
1875  END IF
1876 *
1877 * Do Tests 35 and 36
1878 *
1879  CALL zstt22( n, m, 0, sd, se, d1, dumma, z, ldu, work, m,
1880  $ rwork, result( 35 ) )
1881 *
1882 * Call ZSTEMR to compute D2, do tests.
1883 *
1884 * Compute D2
1885 *
1886  CALL dcopy( n, sd, 1, d5, 1 )
1887  IF( n.GT.0 )
1888  $ CALL dcopy( n-1, se, 1, rwork, 1 )
1889 *
1890  ntest = 37
1891  CALL zstemr( 'N', 'A', n, d5, rwork, vl, vu, il, iu,
1892  $ m, d2, z, ldu, n, iwork( 1 ), tryrac,
1893  $ rwork( n+1 ), lrwork-n, iwork( 2*n+1 ),
1894  $ liwork-2*n, iinfo )
1895  IF( iinfo.NE.0 ) THEN
1896  WRITE( nounit, fmt = 9999 )'ZSTEMR(N,A)', iinfo, n,
1897  $ jtype, ioldsd
1898  info = abs( iinfo )
1899  IF( iinfo.LT.0 ) THEN
1900  RETURN
1901  ELSE
1902  result( 37 ) = ulpinv
1903  GO TO 280
1904  END IF
1905  END IF
1906 *
1907 * Do Test 34
1908 *
1909  temp1 = zero
1910  temp2 = zero
1911 *
1912  DO 260 j = 1, n
1913  temp1 = max( temp1, abs( d1( j ) ), abs( d2( j ) ) )
1914  temp2 = max( temp2, abs( d1( j )-d2( j ) ) )
1915  260 CONTINUE
1916 *
1917  result( 37 ) = temp2 / max( unfl,
1918  $ ulp*max( temp1, temp2 ) )
1919  END IF
1920  270 CONTINUE
1921  280 CONTINUE
1922  ntestt = ntestt + ntest
1923 *
1924 * End of Loop -- Check for RESULT(j) > THRESH
1925 *
1926 *
1927 * Print out tests which fail.
1928 *
1929  DO 290 jr = 1, ntest
1930  IF( result( jr ).GE.thresh ) THEN
1931 *
1932 * If this is the first test to fail,
1933 * print a header to the data file.
1934 *
1935  IF( nerrs.EQ.0 ) THEN
1936  WRITE( nounit, fmt = 9998 )'ZST'
1937  WRITE( nounit, fmt = 9997 )
1938  WRITE( nounit, fmt = 9996 )
1939  WRITE( nounit, fmt = 9995 )'Hermitian'
1940  WRITE( nounit, fmt = 9994 )
1941 *
1942 * Tests performed
1943 *
1944  WRITE( nounit, fmt = 9987 )
1945  END IF
1946  nerrs = nerrs + 1
1947  IF( result( jr ).LT.10000.0d0 ) THEN
1948  WRITE( nounit, fmt = 9989 )n, jtype, ioldsd, jr,
1949  $ result( jr )
1950  ELSE
1951  WRITE( nounit, fmt = 9988 )n, jtype, ioldsd, jr,
1952  $ result( jr )
1953  END IF
1954  END IF
1955  290 CONTINUE
1956  300 CONTINUE
1957  310 CONTINUE
1958 *
1959 * Summary
1960 *
1961  CALL dlasum( 'ZST', nounit, nerrs, ntestt )
1962  RETURN
1963 *
1964  9999 FORMAT( ' ZCHKST: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1965  $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1966 *
1967  9998 FORMAT( / 1x, a3, ' -- Complex Hermitian eigenvalue problem' )
1968  9997 FORMAT( ' Matrix types (see ZCHKST for details): ' )
1969 *
1970  9996 FORMAT( / ' Special Matrices:',
1971  $ / ' 1=Zero matrix. ',
1972  $ ' 5=Diagonal: clustered entries.',
1973  $ / ' 2=Identity matrix. ',
1974  $ ' 6=Diagonal: large, evenly spaced.',
1975  $ / ' 3=Diagonal: evenly spaced entries. ',
1976  $ ' 7=Diagonal: small, evenly spaced.',
1977  $ / ' 4=Diagonal: geometr. spaced entries.' )
1978  9995 FORMAT( ' Dense ', a, ' Matrices:',
1979  $ / ' 8=Evenly spaced eigenvals. ',
1980  $ ' 12=Small, evenly spaced eigenvals.',
1981  $ / ' 9=Geometrically spaced eigenvals. ',
1982  $ ' 13=Matrix with random O(1) entries.',
1983  $ / ' 10=Clustered eigenvalues. ',
1984  $ ' 14=Matrix with large random entries.',
1985  $ / ' 11=Large, evenly spaced eigenvals. ',
1986  $ ' 15=Matrix with small random entries.' )
1987  9994 FORMAT( ' 16=Positive definite, evenly spaced eigenvalues',
1988  $ / ' 17=Positive definite, geometrically spaced eigenvlaues',
1989  $ / ' 18=Positive definite, clustered eigenvalues',
1990  $ / ' 19=Positive definite, small evenly spaced eigenvalues',
1991  $ / ' 20=Positive definite, large evenly spaced eigenvalues',
1992  $ / ' 21=Diagonally dominant tridiagonal, geometrically',
1993  $ ' spaced eigenvalues' )
1994 *
1995  9989 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
1996  $ 4( i4, ',' ), ' result ', i3, ' is', 0p, f8.2 )
1997  9988 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
1998  $ 4( i4, ',' ), ' result ', i3, ' is', 1p, d10.3 )
1999 *
2000  9987 FORMAT( / 'Test performed: see ZCHKST for details.', / )
2001 * End of ZCHKST
2002 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zpteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZPTEQR
Definition: zpteqr.f:147
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine zlatmr(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)
ZLATMR
Definition: zlatmr.f:492
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine zstt22(N, M, KBAND, AD, AE, SD, SE, U, LDU, WORK, LDWORK, RWORK, RESULT)
ZSTT22
Definition: zstt22.f:147
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
double precision function dsxt1(IJOB, D1, N1, D2, N2, ABSTOL, ULP, UNFL)
DSXT1
Definition: dsxt1.f:108
subroutine zhetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
ZHETRD
Definition: zhetrd.f:194
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
subroutine zstt21(N, KBAND, AD, AE, SD, SE, U, LDU, WORK, RWORK, RESULT)
ZSTT21
Definition: zstt21.f:134
subroutine zhpt21(ITYPE, UPLO, N, KBAND, AP, D, E, U, LDU, VP, TAU, WORK, RWORK, RESULT)
ZHPT21
Definition: zhpt21.f:225
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
Definition: zsteqr.f:134
double precision function dlarnd(IDIST, ISEED)
DLARND
Definition: dlarnd.f:75
subroutine zungtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGTR
Definition: zungtr.f:125
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
Definition: zstein.f:184
subroutine zstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZSTEDC
Definition: zstedc.f:215
subroutine zhptrd(UPLO, N, AP, D, E, TAU, INFO)
ZHPTRD
Definition: zhptrd.f:153
subroutine dstech(N, A, B, EIG, TOL, WORK, INFO)
DSTECH
Definition: dstech.f:103
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:45
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zhet21(ITYPE, UPLO, N, KBAND, A, LDA, D, E, U, LDU, V, LDV, TAU, WORK, RWORK, RESULT)
ZHET21
Definition: zhet21.f:213
subroutine zstemr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)
ZSTEMR
Definition: zstemr.f:340
subroutine zupgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
ZUPGTR
Definition: zupgtr.f:116

Here is the call graph for this function:

Here is the caller graph for this function: