LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zdrves()

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.

Definition at line 375 of file zdrves.f.

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