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

◆ cdrves()

subroutine cdrves ( integer  nsizes,
integer, dimension( * )  nn,
integer  ntypes,
logical, dimension( * )  dotype,
integer, dimension( 4 )  iseed,
real  thresh,
integer  nounit,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( lda, * )  h,
complex, dimension( lda, * )  ht,
complex, dimension( * )  w,
complex, dimension( * )  wt,
complex, dimension( ldvs, * )  vs,
integer  ldvs,
real, dimension( 13 )  result,
complex, dimension( * )  work,
integer  nwork,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
logical, dimension( * )  bwork,
integer  info 
)

CDRVES

Purpose:
    CDRVES checks the nonsymmetric eigenvalue (Schur form) problem
    driver CGEES.

    When CDRVES 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,
          CDRVES 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, CDRVES
          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 CDRVES to continue the same random number
          sequence.
[in]THRESH
          THRESH is REAL
          A test will count as "failed" if the "error", computed as
          described above, exceeds THRESH.  Note that the error
          is scaled to be O(1), so THRESH should be a reasonably
          small multiple of 1, e.g., 10 or 100.  In particular,
          it should not depend on the precision (single vs. double)
          or the size of the matrix.  It must be at least zero.
[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 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 array, dimension (LDA, max(NN))
          Another copy of the test matrix A, modified by CGEES.
[out]HT
          HT is COMPLEX array, dimension (LDA, max(NN))
          Yet another copy of the test matrix A, modified by CGEES.
[out]W
          W is COMPLEX array, dimension (max(NN))
          The computed eigenvalues of A.
[out]WT
          WT is COMPLEX array, dimension (max(NN))
          Like W, this array contains the eigenvalues of A,
          but those computed when CGEES only computes a partial
          eigendecomposition, i.e. not Schur vectors
[out]VS
          VS is COMPLEX 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 REAL 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 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 REAL 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  CLATMR, CLATMS, CLATME or CGEES 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 cdrves.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 REAL THRESH
386* ..
387* .. Array Arguments ..
388 LOGICAL BWORK( * ), DOTYPE( * )
389 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
390 REAL RESULT( 13 ), RWORK( * )
391 COMPLEX A( LDA, * ), H( LDA, * ), HT( LDA, * ),
392 $ VS( LDVS, * ), W( * ), WORK( * ), WT( * )
393* ..
394*
395* =====================================================================
396*
397* .. Parameters ..
398 COMPLEX CZERO
399 parameter( czero = ( 0.0e+0, 0.0e+0 ) )
400 COMPLEX CONE
401 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
402 REAL ZERO, ONE
403 parameter( zero = 0.0e+0, one = 1.0e+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,
413 $ NERRS, NFAIL, NMAX, NNWORK, NTEST, NTESTF,
414 $ NTESTT, RSUB, SDIM
415 REAL 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 REAL RES( 2 )
423* ..
424* .. Arrays in Common ..
425 LOGICAL SELVAL( 20 )
426 REAL 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 CSLECT
436 REAL SLAMCH
437 EXTERNAL cslect, slamch
438* ..
439* .. External Subroutines ..
440 EXTERNAL cgees, chst01, clacpy, clatme, clatmr, clatms,
442* ..
443* .. Intrinsic Functions ..
444 INTRINSIC abs, cmplx, 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 ) = 'Complex 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( 'CDRVES', -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 = slamch( 'Safe minimum' )
509 ovfl = one / unfl
510 ulp = slamch( 'Precision' )
511 ulpinv = one / ulp
512 rtulp = sqrt( ulp )
513 rtulpi = one / rtulp
514*
515* Loop over sizes, types
516*
517 nerrs = 0
518*
519 DO 240 jsize = 1, nsizes
520 n = nn( jsize )
521 IF( nsizes.NE.1 ) THEN
522 mtypes = min( maxtyp, ntypes )
523 ELSE
524 mtypes = min( maxtyp+1, ntypes )
525 END IF
526*
527 DO 230 jtype = 1, mtypes
528 IF( .NOT.dotype( jtype ) )
529 $ GO TO 230
530*
531* Save ISEED in case of an error.
532*
533 DO 20 j = 1, 4
534 ioldsd( j ) = iseed( j )
535 20 CONTINUE
536*
537* Compute "A"
538*
539* Control parameters:
540*
541* KMAGN KCONDS KMODE KTYPE
542* =1 O(1) 1 clustered 1 zero
543* =2 large large clustered 2 identity
544* =3 small exponential Jordan
545* =4 arithmetic diagonal, (w/ eigenvalues)
546* =5 random log symmetric, w/ eigenvalues
547* =6 random general, w/ eigenvalues
548* =7 random diagonal
549* =8 random symmetric
550* =9 random general
551* =10 random triangular
552*
553 IF( mtypes.GT.maxtyp )
554 $ GO TO 90
555*
556 itype = ktype( jtype )
557 imode = kmode( jtype )
558*
559* Compute norm
560*
561 GO TO ( 30, 40, 50 )kmagn( jtype )
562*
563 30 CONTINUE
564 anorm = one
565 GO TO 60
566*
567 40 CONTINUE
568 anorm = ovfl*ulp
569 GO TO 60
570*
571 50 CONTINUE
572 anorm = unfl*ulpinv
573 GO TO 60
574*
575 60 CONTINUE
576*
577 CALL claset( 'Full', lda, n, czero, czero, a, lda )
578 iinfo = 0
579 cond = ulpinv
580*
581* Special Matrices -- Identity & Jordan block
582*
583 IF( itype.EQ.1 ) THEN
584*
585* Zero
586*
587 iinfo = 0
588*
589 ELSE IF( itype.EQ.2 ) THEN
590*
591* Identity
592*
593 DO 70 jcol = 1, n
594 a( jcol, jcol ) = cmplx( anorm )
595 70 CONTINUE
596*
597 ELSE IF( itype.EQ.3 ) THEN
598*
599* Jordan Block
600*
601 DO 80 jcol = 1, n
602 a( jcol, jcol ) = cmplx( anorm )
603 IF( jcol.GT.1 )
604 $ a( jcol, jcol-1 ) = cone
605 80 CONTINUE
606*
607 ELSE IF( itype.EQ.4 ) THEN
608*
609* Diagonal Matrix, [Eigen]values Specified
610*
611 CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
612 $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
613 $ iinfo )
614*
615 ELSE IF( itype.EQ.5 ) THEN
616*
617* Symmetric, eigenvalues specified
618*
619 CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
620 $ anorm, n, n, 'N', a, lda, work( n+1 ),
621 $ iinfo )
622*
623 ELSE IF( itype.EQ.6 ) THEN
624*
625* General, eigenvalues specified
626*
627 IF( kconds( jtype ).EQ.1 ) THEN
628 conds = one
629 ELSE IF( kconds( jtype ).EQ.2 ) THEN
630 conds = rtulpi
631 ELSE
632 conds = zero
633 END IF
634*
635 CALL clatme( n, 'D', iseed, work, imode, cond, cone,
636 $ 'T', 'T', 'T', rwork, 4, conds, n, n, anorm,
637 $ a, lda, work( 2*n+1 ), iinfo )
638*
639 ELSE IF( itype.EQ.7 ) THEN
640*
641* Diagonal, random eigenvalues
642*
643 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
644 $ 'T', 'N', work( n+1 ), 1, one,
645 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
646 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
647*
648 ELSE IF( itype.EQ.8 ) THEN
649*
650* Symmetric, random eigenvalues
651*
652 CALL clatmr( n, n, 'D', iseed, 'H', work, 6, one, cone,
653 $ 'T', 'N', work( n+1 ), 1, one,
654 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
655 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
656*
657 ELSE IF( itype.EQ.9 ) THEN
658*
659* General, random eigenvalues
660*
661 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
662 $ 'T', 'N', work( n+1 ), 1, one,
663 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
664 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
665 IF( n.GE.4 ) THEN
666 CALL claset( 'Full', 2, n, czero, czero, a, lda )
667 CALL claset( 'Full', n-3, 1, czero, czero, a( 3, 1 ),
668 $ lda )
669 CALL claset( 'Full', n-3, 2, czero, czero,
670 $ a( 3, n-1 ), lda )
671 CALL claset( 'Full', 1, n, czero, czero, a( n, 1 ),
672 $ lda )
673 END IF
674*
675 ELSE IF( itype.EQ.10 ) THEN
676*
677* Triangular, random eigenvalues
678*
679 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
680 $ 'T', 'N', work( n+1 ), 1, one,
681 $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
682 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
683*
684 ELSE
685*
686 iinfo = 1
687 END IF
688*
689 IF( iinfo.NE.0 ) THEN
690 WRITE( nounit, fmt = 9992 )'Generator', iinfo, n, jtype,
691 $ ioldsd
692 info = abs( iinfo )
693 RETURN
694 END IF
695*
696 90 CONTINUE
697*
698* Test for minimal and generous workspace
699*
700 DO 220 iwk = 1, 2
701 IF( iwk.EQ.1 ) THEN
702 nnwork = 3*n
703 ELSE
704 nnwork = 5*n + 2*n**2
705 END IF
706 nnwork = max( nnwork, 1 )
707*
708* Initialize RESULT
709*
710 DO 100 j = 1, 13
711 result( j ) = -one
712 100 CONTINUE
713*
714* Test with and without sorting of eigenvalues
715*
716 DO 180 isort = 0, 1
717 IF( isort.EQ.0 ) THEN
718 sort = 'N'
719 rsub = 0
720 ELSE
721 sort = 'S'
722 rsub = 6
723 END IF
724*
725* Compute Schur form and Schur vectors, and test them
726*
727 CALL clacpy( 'F', n, n, a, lda, h, lda )
728 CALL cgees( 'V', sort, cslect, n, h, lda, sdim, w, vs,
729 $ ldvs, work, nnwork, rwork, bwork, iinfo )
730 IF( iinfo.NE.0 ) THEN
731 result( 1+rsub ) = ulpinv
732 WRITE( nounit, fmt = 9992 )'CGEES1', iinfo, n,
733 $ jtype, ioldsd
734 info = abs( iinfo )
735 GO TO 190
736 END IF
737*
738* Do Test (1) or Test (7)
739*
740 result( 1+rsub ) = zero
741 DO 120 j = 1, n - 1
742 DO 110 i = j + 1, n
743 IF( h( i, j ).NE.zero )
744 $ result( 1+rsub ) = ulpinv
745 110 CONTINUE
746 120 CONTINUE
747*
748* Do Tests (2) and (3) or Tests (8) and (9)
749*
750 lwork = max( 1, 2*n*n )
751 CALL chst01( n, 1, n, a, lda, h, lda, vs, ldvs, work,
752 $ lwork, rwork, res )
753 result( 2+rsub ) = res( 1 )
754 result( 3+rsub ) = res( 2 )
755*
756* Do Test (4) or Test (10)
757*
758 result( 4+rsub ) = zero
759 DO 130 i = 1, n
760 IF( h( i, i ).NE.w( i ) )
761 $ result( 4+rsub ) = ulpinv
762 130 CONTINUE
763*
764* Do Test (5) or Test (11)
765*
766 CALL clacpy( 'F', n, n, a, lda, ht, lda )
767 CALL cgees( 'N', sort, cslect, n, ht, lda, sdim, wt,
768 $ vs, ldvs, work, nnwork, rwork, bwork,
769 $ iinfo )
770 IF( iinfo.NE.0 ) THEN
771 result( 5+rsub ) = ulpinv
772 WRITE( nounit, fmt = 9992 )'CGEES2', iinfo, n,
773 $ jtype, ioldsd
774 info = abs( iinfo )
775 GO TO 190
776 END IF
777*
778 result( 5+rsub ) = zero
779 DO 150 j = 1, n
780 DO 140 i = 1, n
781 IF( h( i, j ).NE.ht( i, j ) )
782 $ result( 5+rsub ) = ulpinv
783 140 CONTINUE
784 150 CONTINUE
785*
786* Do Test (6) or Test (12)
787*
788 result( 6+rsub ) = zero
789 DO 160 i = 1, n
790 IF( w( i ).NE.wt( i ) )
791 $ result( 6+rsub ) = ulpinv
792 160 CONTINUE
793*
794* Do Test (13)
795*
796 IF( isort.EQ.1 ) THEN
797 result( 13 ) = zero
798 knteig = 0
799 DO 170 i = 1, n
800 IF( cslect( w( i ) ) )
801 $ knteig = knteig + 1
802 IF( i.LT.n ) THEN
803 IF( cslect( w( i+1 ) ) .AND.
804 $ ( .NOT.cslect( w( i ) ) ) )result( 13 )
805 $ = ulpinv
806 END IF
807 170 CONTINUE
808 IF( sdim.NE.knteig )
809 $ result( 13 ) = ulpinv
810 END IF
811*
812 180 CONTINUE
813*
814* End of Loop -- Check for RESULT(j) > THRESH
815*
816 190 CONTINUE
817*
818 ntest = 0
819 nfail = 0
820 DO 200 j = 1, 13
821 IF( result( j ).GE.zero )
822 $ ntest = ntest + 1
823 IF( result( j ).GE.thresh )
824 $ nfail = nfail + 1
825 200 CONTINUE
826*
827 IF( nfail.GT.0 )
828 $ ntestf = ntestf + 1
829 IF( ntestf.EQ.1 ) THEN
830 WRITE( nounit, fmt = 9999 )path
831 WRITE( nounit, fmt = 9998 )
832 WRITE( nounit, fmt = 9997 )
833 WRITE( nounit, fmt = 9996 )
834 WRITE( nounit, fmt = 9995 )thresh
835 WRITE( nounit, fmt = 9994 )
836 ntestf = 2
837 END IF
838*
839 DO 210 j = 1, 13
840 IF( result( j ).GE.thresh ) THEN
841 WRITE( nounit, fmt = 9993 )n, iwk, ioldsd, jtype,
842 $ j, result( j )
843 END IF
844 210 CONTINUE
845*
846 nerrs = nerrs + nfail
847 ntestt = ntestt + ntest
848*
849 220 CONTINUE
850 230 CONTINUE
851 240 CONTINUE
852*
853* Summary
854*
855 CALL slasum( path, nounit, nerrs, ntestt )
856*
857 9999 FORMAT( / 1x, a3, ' -- Complex Schur Form Decomposition Driver',
858 $ / ' Matrix types (see CDRVES for details): ' )
859*
860 9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
861 $ ' ', ' 5=Diagonal: geometr. spaced entries.',
862 $ / ' 2=Identity matrix. ', ' 6=Diagona',
863 $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
864 $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
865 $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
866 $ 'mall, evenly spaced.' )
867 9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
868 $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
869 $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
870 $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
871 $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
872 $ 'lex ', a6, / ' 12=Well-cond., random complex ', a6, ' ',
873 $ ' 17=Ill-cond., large rand. complx ', a4, / ' 13=Ill-condi',
874 $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
875 $ ' complx ', a4 )
876 9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
877 $ 'with small random entries.', / ' 20=Matrix with large ran',
878 $ 'dom entries. ', / )
879 9995 FORMAT( ' Tests performed with test threshold =', f8.2,
880 $ / ' ( A denotes A on input and T denotes A on output)',
881 $ / / ' 1 = 0 if T in Schur form (no sort), ',
882 $ ' 1/ulp otherwise', /
883 $ ' 2 = | A - VS T transpose(VS) | / ( n |A| ulp ) (no sort)',
884 $ / ' 3 = | I - VS transpose(VS) | / ( n ulp ) (no sort) ',
885 $ / ' 4 = 0 if W are eigenvalues of T (no sort),',
886 $ ' 1/ulp otherwise', /
887 $ ' 5 = 0 if T same no matter if VS computed (no sort),',
888 $ ' 1/ulp otherwise', /
889 $ ' 6 = 0 if W same no matter if VS computed (no sort)',
890 $ ', 1/ulp otherwise' )
891 9994 FORMAT( ' 7 = 0 if T in Schur form (sort), ', ' 1/ulp otherwise',
892 $ / ' 8 = | A - VS T transpose(VS) | / ( n |A| ulp ) (sort)',
893 $ / ' 9 = | I - VS transpose(VS) | / ( n ulp ) (sort) ',
894 $ / ' 10 = 0 if W are eigenvalues of T (sort),',
895 $ ' 1/ulp otherwise', /
896 $ ' 11 = 0 if T same no matter if VS computed (sort),',
897 $ ' 1/ulp otherwise', /
898 $ ' 12 = 0 if W same no matter if VS computed (sort),',
899 $ ' 1/ulp otherwise', /
900 $ ' 13 = 0 if sorting successful, 1/ulp otherwise', / )
901 9993 FORMAT( ' N=', i5, ', IWK=', i2, ', seed=', 4( i4, ',' ),
902 $ ' type ', i2, ', test(', i2, ')=', g10.3 )
903 9992 FORMAT( ' CDRVES: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
904 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
905*
906 RETURN
907*
908* End of CDRVES
909*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chst01(n, ilo, ihi, a, lda, h, ldh, q, ldq, work, lwork, rwork, result)
CHST01
Definition chst01.f:140
subroutine clatme(n, dist, iseed, d, mode, cond, dmax, rsign, upper, sim, ds, modes, conds, kl, ku, anorm, a, lda, work, info)
CLATME
Definition clatme.f:301
subroutine clatmr(m, n, dist, iseed, sym, d, mode, cond, dmax, rsign, grade, dl, model, condl, dr, moder, condr, pivtng, ipivot, kl, ku, sparse, anorm, pack, a, lda, iwork, info)
CLATMR
Definition clatmr.f:490
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
logical function cslect(z)
CSLECT
Definition cslect.f:56
subroutine cgees(jobvs, sort, select, n, a, lda, sdim, w, vs, ldvs, work, lwork, rwork, bwork, info)
CGEES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE m...
Definition cgees.f:197
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:106
subroutine slasum(type, iounit, ie, nrun)
SLASUM
Definition slasum.f:41
Here is the call graph for this function:
Here is the caller graph for this function: