LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine dchkst2stg ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
integer  NOUNIT,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, 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,
double precision, dimension( ldu, * )  U,
integer  LDU,
double precision, dimension( ldu, * )  V,
double precision, dimension( * )  VP,
double precision, dimension( * )  TAU,
double precision, dimension( ldu, * )  Z,
double precision, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

DCHKST2STG

Purpose:
 DCHKST2STG  checks the symmetric eigenvalue problem routines
 using the 2-stage reduction techniques. Since the generation
 of Q or the vectors is not available in this release, we only 
 compare the eigenvalue resulting when using the 2-stage to the 
 one considered as reference using the standard 1-stage reduction
 DSYTRD. For that, we call the standard DSYTRD and compute D1 using 
 DSTEQR, then we call the 2-stage DSYTRD_2STAGE with Upper and Lower
 and we compute D2 and D3 using DSTEQR and then we replaced tests
 3 and 4 by tests 11 and 12. test 1 and 2 remain to verify that 
 the 1-stage results are OK and can be trusted.
 This testing routine will converge to the DCHKST in the next 
 release when vectors and generation of Q will be implemented.

    DSYTRD factors A as  U S U' , where ' means transpose,
    S is symmetric tridiagonal, and U is orthogonal.
    DSYTRD can use either just the lower or just the upper triangle
    of A; DCHKST2STG 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.

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

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

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

    DSTEQR factors S as  Z D1 Z' , where Z is the orthogonal
    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.

    DPTEQR factors S as  Z4 D4 Z4' , for a
    symmetric 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.

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

    DSTEDC factors S as Z D1 Z' , where Z is the orthogonal
    matrix of eigenvectors and D1 is a diagonal matrix with
    the eigenvalues on the diagonal ('I' option). It may also
    update an input orthogonal matrix, usually the output
    from DSYTRD/DORGTR or DSPTRD/DOPGTR ('V' option). It may
    also just compute eigenvalues ('N' option).

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

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

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

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

 (3)     | A - V S V' | / ( |A| n ulp ) DSYTRD( UPLO='L', ... )
         replaced by | D1 - D2 | / ( |D1| ulp ) where D1 is the 
         eigenvalue matrix computed using S and D2 is the 
         eigenvalue matrix computed using S_2stage the output of
         DSYTRD_2STAGE("N", "U",....). D1 and D2 are computed 
         via DSTEQR('N',...)  

 (4)     | I - UV' | / ( n ulp )        DORGTR( UPLO='L', ... )
         replaced by | D1 - D3 | / ( |D1| ulp ) where D1 is the 
         eigenvalue matrix computed using S and D3 is the 
         eigenvalue matrix computed using S_2stage the output of
         DSYTRD_2STAGE("N", "L",....). D1 and D3 are computed 
         via DSTEQR('N',...)  

 (5-8)   Same as 1-4, but for DSPTRD and DOPGTR.

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

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

 (11)    | D1 - D2 | / ( |D1| ulp )        DSTEQR('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 ) DPTEQR('V',...)

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

 (16)    | D4 - D5 | / ( 100 |D4| ulp )       DPTEQR('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, SSTEIN

 (21)    | I - Y Y' | / ( n ulp )          DSTEBZ, SSTEIN

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

 (23)    | I - ZZ' | / ( n ulp )           DSTEDC('I')

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

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

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

 Test 27 is disabled at the moment because DSTEMR 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
                                              DSTEMR('V', 'A')

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

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

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

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

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

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

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

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

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

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

 (37)    ( max { min | WA2(i)-WA3(j) | } +
            i     j
           max { min | WA3(i)-WA2(j) | } ) / ( |D3| ulp )
            i     j
         DSTEMR('N', 'A') vs. SSTEMR('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 orthogonal 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 orthogonal 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 orthogonal 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) 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,
          DCHKST2STG 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, DCHKST2STG
          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 DCHKST2STG 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DSYTRD.
          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
          DSYTRD.  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 DSTEQR 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 DSTEQR 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 DPTEQR(V).
          DPTEQR 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 DPTEQR(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 DOUBLE PRECISION array of
                             dimension( LDU, max(NN) ).
          The orthogonal matrix computed by DSYTRD + DORGTR.
[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 DOUBLE PRECISION array of
                             dimension( LDU, max(NN) ).
          The Housholder vectors computed by DSYTRD 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 DSYTRD, 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 DORGTR, set those entries to
          1 before using them, and then restore them later.)
[out]VP
          VP is DOUBLE PRECISION array of
                      dimension( max(NN)*max(NN+1)/2 )
          The matrix V stored in packed format.
[out]TAU
          TAU is DOUBLE PRECISION array of
                             dimension( max(NN) )
          The Householder factors computed by DSYTRD in reducing A
          to tridiagonal form.
[out]Z
          Z is DOUBLE PRECISION array of
                             dimension( LDU, max(NN) ).
          The orthogonal matrix of eigenvectors computed by DSTEQR,
          DPTEQR, and DSTEIN.
[out]WORK
          WORK is DOUBLE PRECISION 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]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  DLATMR, SLATMS, DSYTRD, DORGTR, DSTEQR, SSTERF,
              or DORMC2 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
December 2016

Definition at line 614 of file dchkst2stg.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: