LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrves ( 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( lda, * )  H,
complex*16, dimension( lda, * )  HT,
complex*16, dimension( * )  W,
complex*16, dimension( * )  WT,
complex*16, dimension( ldvs, * )  VS,
integer  LDVS,
double precision, dimension( 13 )  RESULT,
complex*16, dimension( * )  WORK,
integer  NWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
logical, dimension( * )  BWORK,
integer  INFO 
)

ZDRVES

Purpose:
    ZDRVES checks the nonsymmetric eigenvalue (Schur form) problem
    driver ZGEES.

    When ZDRVES 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 nonsymmetric eigenroutines.  For each matrix, 13
    tests will be performed:

    (1)     0 if T is in Schur form, 1/ulp otherwise
           (no sorting of eigenvalues)

    (2)     | A - VS T VS' | / ( n |A| ulp )

      Here VS is the matrix of Schur eigenvectors, and T is in Schur
      form  (no sorting of eigenvalues).

    (3)     | I - VS VS' | / ( n ulp ) (no sorting of eigenvalues).

    (4)     0     if W are eigenvalues of T
            1/ulp otherwise
            (no sorting of eigenvalues)

    (5)     0     if T(with VS) = T(without VS),
            1/ulp otherwise
            (no sorting of eigenvalues)

    (6)     0     if eigenvalues(with VS) = eigenvalues(without VS),
            1/ulp otherwise
            (no sorting of eigenvalues)

    (7)     0 if T is in Schur form, 1/ulp otherwise
            (with sorting of eigenvalues)

    (8)     | A - VS T VS' | / ( n |A| ulp )

      Here VS is the matrix of Schur eigenvectors, and T is in Schur
      form  (with sorting of eigenvalues).

    (9)     | I - VS VS' | / ( n ulp ) (with sorting of eigenvalues).

    (10)    0     if W are eigenvalues of T
            1/ulp otherwise
            (with sorting of eigenvalues)

    (11)    0     if T(with VS) = T(without VS),
            1/ulp otherwise
            (with sorting of eigenvalues)

    (12)    0     if eigenvalues(with VS) = eigenvalues(without VS),
            1/ulp otherwise
            (with sorting of eigenvalues)

    (13)    if sorting worked and SDIM is the number of
            eigenvalues which were SELECTed

    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 (transposed) Jordan block, with 1's on the diagonal.

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

    (7)  Same as (4), but multiplied by a constant near
         the overflow threshold
    (8)  Same as (4), but multiplied by a constant near
         the underflow threshold

    (9)  A matrix of the form  U' T U, where U is unitary and
         T has evenly spaced entries 1, ..., ULP with random
         complex angles on the diagonal and random O(1) entries in
         the upper triangle.

    (10) A matrix of the form  U' T U, where U is unitary and
         T has geometrically spaced entries 1, ..., ULP with random
         complex angles on the diagonal and random O(1) entries in
         the upper triangle.

    (11) A matrix of the form  U' T U, where U is orthogonal and
         T has "clustered" entries 1, ULP,..., ULP with random
         complex angles on the diagonal and random O(1) entries in
         the upper triangle.

    (12) A matrix of the form  U' T U, where U is unitary and
         T has complex eigenvalues randomly chosen from
         ULP < |z| < 1   and random O(1) entries in the upper
         triangle.

    (13) A matrix of the form  X' T X, where X has condition
         SQRT( ULP ) and T has evenly spaced entries 1, ..., ULP
         with random complex angles on the diagonal and random O(1)
         entries in the upper triangle.

    (14) A matrix of the form  X' T X, where X has condition
         SQRT( ULP ) and T has geometrically spaced entries
         1, ..., ULP with random complex angles on the diagonal
         and random O(1) entries in the upper triangle.

    (15) A matrix of the form  X' T X, where X has condition
         SQRT( ULP ) and T has "clustered" entries 1, ULP,..., ULP
         with random complex angles on the diagonal and random O(1)
         entries in the upper triangle.

    (16) A matrix of the form  X' T X, where X has condition
         SQRT( ULP ) and T has complex eigenvalues randomly chosen
         from ULP < |z| < 1 and random O(1) entries in the upper
         triangle.

    (17) Same as (16), but multiplied by a constant
         near the overflow threshold
    (18) Same as (16), but multiplied by a constant
         near the underflow threshold

    (19) Nonsymmetric matrix with random entries chosen from (-1,1).
         If N is at least 4, all entries in first two rows and last
         row, and first column and last two columns are zero.
    (20) Same as (19), but multiplied by a constant
         near the overflow threshold
    (21) Same as (19), but multiplied by a constant
         near the underflow threshold
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of sizes of matrices to use.  If it is zero,
          ZDRVES 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, ZDRVES
          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 ZDRVES 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 INFO not equal to 0.)
[out]A
          A is COMPLEX*16 array, dimension (LDA, max(NN))
          Used to hold the matrix whose eigenvalues are to be
          computed.  On exit, A contains the last matrix actually used.
[in]LDA
          LDA is INTEGER
          The leading dimension of A, and H. LDA must be at
          least 1 and at least max( NN ).
[out]H
          H is COMPLEX*16 array, dimension (LDA, max(NN))
          Another copy of the test matrix A, modified by ZGEES.
[out]HT
          HT is COMPLEX*16 array, dimension (LDA, max(NN))
          Yet another copy of the test matrix A, modified by ZGEES.
[out]W
          W is COMPLEX*16 array, dimension (max(NN))
          The computed eigenvalues of A.
[out]WT
          WT is COMPLEX*16 array, dimension (max(NN))
          Like W, this array contains the eigenvalues of A,
          but those computed when ZGEES only computes a partial
          eigendecomposition, i.e. not Schur vectors
[out]VS
          VS is COMPLEX*16 array, dimension (LDVS, max(NN))
          VS holds the computed Schur vectors.
[in]LDVS
          LDVS is INTEGER
          Leading dimension of VS. Must be at least max(1,max(NN)).
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (13)
          The values computed by the 13 tests described above.
          The values are currently limited to 1/ulp, to avoid overflow.
[out]WORK
          WORK is COMPLEX*16 array, dimension (NWORK)
[in]NWORK
          NWORK is INTEGER
          The number of entries in WORK.  This must be at least
          5*NN(j)+2*NN(j)**2 for all j.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(NN))
[out]IWORK
          IWORK is INTEGER array, dimension (max(NN))
[out]BWORK
          BWORK is LOGICAL array, dimension (max(NN))
[out]INFO
          INFO is INTEGER
          If 0, then everything ran OK.
           -1: NSIZES < 0
           -2: Some NN(j) < 0
           -3: NTYPES < 0
           -6: THRESH < 0
           -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
          -15: LDVS < 1 or LDVS < NMAX, where NMAX is max( NN(j) ).
          -18: NWORK too small.
          If  ZLATMR, CLATMS, CLATME or ZGEES 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.
     NMAX            Largest value in NN.
     NERRS           The number of tests which have exceeded THRESH
     COND, CONDS,
     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.
     RTULP, RTULPI   Square roots of the previous 4 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) )
     KCONDS(j)       Select whether CONDS is to be 1 or
                     1/sqrt(ulp).  (0 means irrelevant.)
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
June 2016

Definition at line 380 of file zdrves.f.

380 *
381 * -- LAPACK test routine (version 3.6.1) --
382 * -- LAPACK is a software package provided by Univ. of Tennessee, --
383 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
384 * June 2016
385 *
386 * .. Scalar Arguments ..
387  INTEGER info, lda, ldvs, nounit, nsizes, ntypes, nwork
388  DOUBLE PRECISION thresh
389 * ..
390 * .. Array Arguments ..
391  LOGICAL bwork( * ), dotype( * )
392  INTEGER iseed( 4 ), iwork( * ), nn( * )
393  DOUBLE PRECISION result( 13 ), rwork( * )
394  COMPLEX*16 a( lda, * ), h( lda, * ), ht( lda, * ),
395  $ vs( ldvs, * ), w( * ), work( * ), wt( * )
396 * ..
397 *
398 * =====================================================================
399 *
400 * .. Parameters ..
401  COMPLEX*16 czero
402  parameter ( czero = ( 0.0d+0, 0.0d+0 ) )
403  COMPLEX*16 cone
404  parameter ( cone = ( 1.0d+0, 0.0d+0 ) )
405  DOUBLE PRECISION zero, one
406  parameter ( zero = 0.0d+0, one = 1.0d+0 )
407  INTEGER maxtyp
408  parameter ( maxtyp = 21 )
409 * ..
410 * .. Local Scalars ..
411  LOGICAL badnn
412  CHARACTER sort
413  CHARACTER*3 path
414  INTEGER i, iinfo, imode, isort, itype, iwk, j, jcol,
415  $ jsize, jtype, knteig, lwork, mtypes, n, nerrs,
416  $ nfail, nmax, nnwork, ntest, ntestf, ntestt,
417  $ rsub, sdim
418  DOUBLE PRECISION anorm, cond, conds, ovfl, rtulp, rtulpi, ulp,
419  $ ulpinv, unfl
420 * ..
421 * .. Local Arrays ..
422  INTEGER idumma( 1 ), ioldsd( 4 ), kconds( maxtyp ),
423  $ kmagn( maxtyp ), kmode( maxtyp ),
424  $ ktype( maxtyp )
425  DOUBLE PRECISION res( 2 )
426 * ..
427 * .. Arrays in Common ..
428  LOGICAL selval( 20 )
429  DOUBLE PRECISION selwi( 20 ), selwr( 20 )
430 * ..
431 * .. Scalars in Common ..
432  INTEGER seldim, selopt
433 * ..
434 * .. Common blocks ..
435  COMMON / sslct / selopt, seldim, selval, selwr, selwi
436 * ..
437 * .. External Functions ..
438  LOGICAL zslect
439  DOUBLE PRECISION dlamch
440  EXTERNAL zslect, dlamch
441 * ..
442 * .. External Subroutines ..
443  EXTERNAL dlabad, dlasum, xerbla, zgees, zhst01, zlacpy,
445 * ..
446 * .. Intrinsic Functions ..
447  INTRINSIC abs, dcmplx, max, min, sqrt
448 * ..
449 * .. Data statements ..
450  DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
451  DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
452  $ 3, 1, 2, 3 /
453  DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
454  $ 1, 5, 5, 5, 4, 3, 1 /
455  DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
456 * ..
457 * .. Executable Statements ..
458 *
459  path( 1: 1 ) = 'Zomplex precision'
460  path( 2: 3 ) = 'ES'
461 *
462 * Check for errors
463 *
464  ntestt = 0
465  ntestf = 0
466  info = 0
467  selopt = 0
468 *
469 * Important constants
470 *
471  badnn = .false.
472  nmax = 0
473  DO 10 j = 1, nsizes
474  nmax = max( nmax, nn( j ) )
475  IF( nn( j ).LT.0 )
476  $ badnn = .true.
477  10 CONTINUE
478 *
479 * Check for errors
480 *
481  IF( nsizes.LT.0 ) THEN
482  info = -1
483  ELSE IF( badnn ) THEN
484  info = -2
485  ELSE IF( ntypes.LT.0 ) THEN
486  info = -3
487  ELSE IF( thresh.LT.zero ) THEN
488  info = -6
489  ELSE IF( nounit.LE.0 ) THEN
490  info = -7
491  ELSE IF( lda.LT.1 .OR. lda.LT.nmax ) THEN
492  info = -9
493  ELSE IF( ldvs.LT.1 .OR. ldvs.LT.nmax ) THEN
494  info = -15
495  ELSE IF( 5*nmax+2*nmax**2.GT.nwork ) THEN
496  info = -18
497  END IF
498 *
499  IF( info.NE.0 ) THEN
500  CALL xerbla( 'ZDRVES', -info )
501  RETURN
502  END IF
503 *
504 * Quick return if nothing to do
505 *
506  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
507  $ RETURN
508 *
509 * More Important constants
510 *
511  unfl = dlamch( 'Safe minimum' )
512  ovfl = one / unfl
513  CALL dlabad( unfl, ovfl )
514  ulp = dlamch( 'Precision' )
515  ulpinv = one / ulp
516  rtulp = sqrt( ulp )
517  rtulpi = one / rtulp
518 *
519 * Loop over sizes, types
520 *
521  nerrs = 0
522 *
523  DO 240 jsize = 1, nsizes
524  n = nn( jsize )
525  IF( nsizes.NE.1 ) THEN
526  mtypes = min( maxtyp, ntypes )
527  ELSE
528  mtypes = min( maxtyp+1, ntypes )
529  END IF
530 *
531  DO 230 jtype = 1, mtypes
532  IF( .NOT.dotype( jtype ) )
533  $ GO TO 230
534 *
535 * Save ISEED in case of an error.
536 *
537  DO 20 j = 1, 4
538  ioldsd( j ) = iseed( j )
539  20 CONTINUE
540 *
541 * Compute "A"
542 *
543 * Control parameters:
544 *
545 * KMAGN KCONDS KMODE KTYPE
546 * =1 O(1) 1 clustered 1 zero
547 * =2 large large clustered 2 identity
548 * =3 small exponential Jordan
549 * =4 arithmetic diagonal, (w/ eigenvalues)
550 * =5 random log symmetric, w/ eigenvalues
551 * =6 random general, w/ eigenvalues
552 * =7 random diagonal
553 * =8 random symmetric
554 * =9 random general
555 * =10 random triangular
556 *
557  IF( mtypes.GT.maxtyp )
558  $ GO TO 90
559 *
560  itype = ktype( jtype )
561  imode = kmode( jtype )
562 *
563 * Compute norm
564 *
565  GO TO ( 30, 40, 50 )kmagn( jtype )
566 *
567  30 CONTINUE
568  anorm = one
569  GO TO 60
570 *
571  40 CONTINUE
572  anorm = ovfl*ulp
573  GO TO 60
574 *
575  50 CONTINUE
576  anorm = unfl*ulpinv
577  GO TO 60
578 *
579  60 CONTINUE
580 *
581  CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
582  iinfo = 0
583  cond = ulpinv
584 *
585 * Special Matrices -- Identity & Jordan block
586 *
587  IF( itype.EQ.1 ) THEN
588 *
589 * Zero
590 *
591  iinfo = 0
592 *
593  ELSE IF( itype.EQ.2 ) THEN
594 *
595 * Identity
596 *
597  DO 70 jcol = 1, n
598  a( jcol, jcol ) = dcmplx( anorm )
599  70 CONTINUE
600 *
601  ELSE IF( itype.EQ.3 ) THEN
602 *
603 * Jordan Block
604 *
605  DO 80 jcol = 1, n
606  a( jcol, jcol ) = dcmplx( anorm )
607  IF( jcol.GT.1 )
608  $ a( jcol, jcol-1 ) = cone
609  80 CONTINUE
610 *
611  ELSE IF( itype.EQ.4 ) THEN
612 *
613 * Diagonal Matrix, [Eigen]values Specified
614 *
615  CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
616  $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
617  $ iinfo )
618 *
619  ELSE IF( itype.EQ.5 ) THEN
620 *
621 * Symmetric, eigenvalues specified
622 *
623  CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
624  $ anorm, n, n, 'N', a, lda, work( n+1 ),
625  $ iinfo )
626 *
627  ELSE IF( itype.EQ.6 ) THEN
628 *
629 * General, eigenvalues specified
630 *
631  IF( kconds( jtype ).EQ.1 ) THEN
632  conds = one
633  ELSE IF( kconds( jtype ).EQ.2 ) THEN
634  conds = rtulpi
635  ELSE
636  conds = zero
637  END IF
638 *
639  CALL zlatme( n, 'D', iseed, work, imode, cond, cone,
640  $ 'T', 'T', 'T', rwork, 4, conds, n, n, anorm,
641  $ a, lda, work( 2*n+1 ), iinfo )
642 *
643  ELSE IF( itype.EQ.7 ) THEN
644 *
645 * Diagonal, random eigenvalues
646 *
647  CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
648  $ 'T', 'N', work( n+1 ), 1, one,
649  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
650  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
651 *
652  ELSE IF( itype.EQ.8 ) THEN
653 *
654 * Symmetric, random eigenvalues
655 *
656  CALL zlatmr( n, n, 'D', iseed, 'H', work, 6, one, cone,
657  $ 'T', 'N', work( n+1 ), 1, one,
658  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
659  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
660 *
661  ELSE IF( itype.EQ.9 ) THEN
662 *
663 * General, random eigenvalues
664 *
665  CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
666  $ 'T', 'N', work( n+1 ), 1, one,
667  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
668  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
669  IF( n.GE.4 ) THEN
670  CALL zlaset( 'Full', 2, n, czero, czero, a, lda )
671  CALL zlaset( 'Full', n-3, 1, czero, czero, a( 3, 1 ),
672  $ lda )
673  CALL zlaset( 'Full', n-3, 2, czero, czero,
674  $ a( 3, n-1 ), lda )
675  CALL zlaset( 'Full', 1, n, czero, czero, a( n, 1 ),
676  $ lda )
677  END IF
678 *
679  ELSE IF( itype.EQ.10 ) THEN
680 *
681 * Triangular, random eigenvalues
682 *
683  CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
684  $ 'T', 'N', work( n+1 ), 1, one,
685  $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
686  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
687 *
688  ELSE
689 *
690  iinfo = 1
691  END IF
692 *
693  IF( iinfo.NE.0 ) THEN
694  WRITE( nounit, fmt = 9992 )'Generator', iinfo, n, jtype,
695  $ ioldsd
696  info = abs( iinfo )
697  RETURN
698  END IF
699 *
700  90 CONTINUE
701 *
702 * Test for minimal and generous workspace
703 *
704  DO 220 iwk = 1, 2
705  IF( iwk.EQ.1 ) THEN
706  nnwork = 3*n
707  ELSE
708  nnwork = 5*n + 2*n**2
709  END IF
710  nnwork = max( nnwork, 1 )
711 *
712 * Initialize RESULT
713 *
714  DO 100 j = 1, 13
715  result( j ) = -one
716  100 CONTINUE
717 *
718 * Test with and without sorting of eigenvalues
719 *
720  DO 180 isort = 0, 1
721  IF( isort.EQ.0 ) THEN
722  sort = 'N'
723  rsub = 0
724  ELSE
725  sort = 'S'
726  rsub = 6
727  END IF
728 *
729 * Compute Schur form and Schur vectors, and test them
730 *
731  CALL zlacpy( 'F', n, n, a, lda, h, lda )
732  CALL zgees( 'V', sort, zslect, n, h, lda, sdim, w, vs,
733  $ ldvs, work, nnwork, rwork, bwork, iinfo )
734  IF( iinfo.NE.0 ) THEN
735  result( 1+rsub ) = ulpinv
736  WRITE( nounit, fmt = 9992 )'ZGEES1', iinfo, n,
737  $ jtype, ioldsd
738  info = abs( iinfo )
739  GO TO 190
740  END IF
741 *
742 * Do Test (1) or Test (7)
743 *
744  result( 1+rsub ) = zero
745  DO 120 j = 1, n - 1
746  DO 110 i = j + 1, n
747  IF( h( i, j ).NE.zero )
748  $ result( 1+rsub ) = ulpinv
749  110 CONTINUE
750  120 CONTINUE
751 *
752 * Do Tests (2) and (3) or Tests (8) and (9)
753 *
754  lwork = max( 1, 2*n*n )
755  CALL zhst01( n, 1, n, a, lda, h, lda, vs, ldvs, work,
756  $ lwork, rwork, res )
757  result( 2+rsub ) = res( 1 )
758  result( 3+rsub ) = res( 2 )
759 *
760 * Do Test (4) or Test (10)
761 *
762  result( 4+rsub ) = zero
763  DO 130 i = 1, n
764  IF( h( i, i ).NE.w( i ) )
765  $ result( 4+rsub ) = ulpinv
766  130 CONTINUE
767 *
768 * Do Test (5) or Test (11)
769 *
770  CALL zlacpy( 'F', n, n, a, lda, ht, lda )
771  CALL zgees( 'N', sort, zslect, n, ht, lda, sdim, wt,
772  $ vs, ldvs, work, nnwork, rwork, bwork,
773  $ iinfo )
774  IF( iinfo.NE.0 ) THEN
775  result( 5+rsub ) = ulpinv
776  WRITE( nounit, fmt = 9992 )'ZGEES2', iinfo, n,
777  $ jtype, ioldsd
778  info = abs( iinfo )
779  GO TO 190
780  END IF
781 *
782  result( 5+rsub ) = zero
783  DO 150 j = 1, n
784  DO 140 i = 1, n
785  IF( h( i, j ).NE.ht( i, j ) )
786  $ result( 5+rsub ) = ulpinv
787  140 CONTINUE
788  150 CONTINUE
789 *
790 * Do Test (6) or Test (12)
791 *
792  result( 6+rsub ) = zero
793  DO 160 i = 1, n
794  IF( w( i ).NE.wt( i ) )
795  $ result( 6+rsub ) = ulpinv
796  160 CONTINUE
797 *
798 * Do Test (13)
799 *
800  IF( isort.EQ.1 ) THEN
801  result( 13 ) = zero
802  knteig = 0
803  DO 170 i = 1, n
804  IF( zslect( w( i ) ) )
805  $ knteig = knteig + 1
806  IF( i.LT.n ) THEN
807  IF( zslect( w( i+1 ) ) .AND.
808  $ ( .NOT.zslect( w( i ) ) ) )result( 13 )
809  $ = ulpinv
810  END IF
811  170 CONTINUE
812  IF( sdim.NE.knteig )
813  $ result( 13 ) = ulpinv
814  END IF
815 *
816  180 CONTINUE
817 *
818 * End of Loop -- Check for RESULT(j) > THRESH
819 *
820  190 CONTINUE
821 *
822  ntest = 0
823  nfail = 0
824  DO 200 j = 1, 13
825  IF( result( j ).GE.zero )
826  $ ntest = ntest + 1
827  IF( result( j ).GE.thresh )
828  $ nfail = nfail + 1
829  200 CONTINUE
830 *
831  IF( nfail.GT.0 )
832  $ ntestf = ntestf + 1
833  IF( ntestf.EQ.1 ) THEN
834  WRITE( nounit, fmt = 9999 )path
835  WRITE( nounit, fmt = 9998 )
836  WRITE( nounit, fmt = 9997 )
837  WRITE( nounit, fmt = 9996 )
838  WRITE( nounit, fmt = 9995 )thresh
839  WRITE( nounit, fmt = 9994 )
840  ntestf = 2
841  END IF
842 *
843  DO 210 j = 1, 13
844  IF( result( j ).GE.thresh ) THEN
845  WRITE( nounit, fmt = 9993 )n, iwk, ioldsd, jtype,
846  $ j, result( j )
847  END IF
848  210 CONTINUE
849 *
850  nerrs = nerrs + nfail
851  ntestt = ntestt + ntest
852 *
853  220 CONTINUE
854  230 CONTINUE
855  240 CONTINUE
856 *
857 * Summary
858 *
859  CALL dlasum( path, nounit, nerrs, ntestt )
860 *
861  9999 FORMAT( / 1x, a3, ' -- Complex Schur Form Decomposition Driver',
862  $ / ' Matrix types (see ZDRVES for details): ' )
863 *
864  9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
865  $ ' ', ' 5=Diagonal: geometr. spaced entries.',
866  $ / ' 2=Identity matrix. ', ' 6=Diagona',
867  $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
868  $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
869  $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
870  $ 'mall, evenly spaced.' )
871  9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
872  $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
873  $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
874  $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
875  $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
876  $ 'lex ', a6, / ' 12=Well-cond., random complex ', a6, ' ',
877  $ ' 17=Ill-cond., large rand. complx ', a4, / ' 13=Ill-condi',
878  $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
879  $ ' complx ', a4 )
880  9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
881  $ 'with small random entries.', / ' 20=Matrix with large ran',
882  $ 'dom entries. ', / )
883  9995 FORMAT( ' Tests performed with test threshold =', f8.2,
884  $ / ' ( A denotes A on input and T denotes A on output)',
885  $ / / ' 1 = 0 if T in Schur form (no sort), ',
886  $ ' 1/ulp otherwise', /
887  $ ' 2 = | A - VS T transpose(VS) | / ( n |A| ulp ) (no sort)',
888  $ / ' 3 = | I - VS transpose(VS) | / ( n ulp ) (no sort) ',
889  $ / ' 4 = 0 if W are eigenvalues of T (no sort),',
890  $ ' 1/ulp otherwise', /
891  $ ' 5 = 0 if T same no matter if VS computed (no sort),',
892  $ ' 1/ulp otherwise', /
893  $ ' 6 = 0 if W same no matter if VS computed (no sort)',
894  $ ', 1/ulp otherwise' )
895  9994 FORMAT( ' 7 = 0 if T in Schur form (sort), ', ' 1/ulp otherwise',
896  $ / ' 8 = | A - VS T transpose(VS) | / ( n |A| ulp ) (sort)',
897  $ / ' 9 = | I - VS transpose(VS) | / ( n ulp ) (sort) ',
898  $ / ' 10 = 0 if W are eigenvalues of T (sort),',
899  $ ' 1/ulp otherwise', /
900  $ ' 11 = 0 if T same no matter if VS computed (sort),',
901  $ ' 1/ulp otherwise', /
902  $ ' 12 = 0 if W same no matter if VS computed (sort),',
903  $ ' 1/ulp otherwise', /
904  $ ' 13 = 0 if sorting successful, 1/ulp otherwise', / )
905  9993 FORMAT( ' N=', i5, ', IWK=', i2, ', seed=', 4( i4, ',' ),
906  $ ' type ', i2, ', test(', i2, ')=', g10.3 )
907  9992 FORMAT( ' ZDRVES: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
908  $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
909 *
910  RETURN
911 *
912 * End of ZDRVES
913 *
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 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
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgees(JOBVS, SORT, SELECT, N, A, LDA, SDIM, W, VS, LDVS, WORK, LWORK, RWORK, BWORK, INFO)
ZGEES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition: zgees.f:199
subroutine zhst01(N, ILO, IHI, A, LDA, H, LDH, Q, LDQ, WORK, LWORK, RWORK, RESULT)
ZHST01
Definition: zhst01.f:142
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 xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zlatme(N, DIST, ISEED, D, MODE, COND, DMAX, RSIGN, UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM, A, LDA, WORK, INFO)
ZLATME
Definition: zlatme.f:303
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:45
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
logical function zslect(Z)
ZSLECT
Definition: zslect.f:58

Here is the call graph for this function:

Here is the caller graph for this function: