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

◆ zdrvev()

subroutine zdrvev ( 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( * )  W,
complex*16, dimension( * )  W1,
complex*16, dimension( ldvl, * )  VL,
integer  LDVL,
complex*16, dimension( ldvr, * )  VR,
integer  LDVR,
complex*16, dimension( ldlre, * )  LRE,
integer  LDLRE,
double precision, dimension( 7 )  RESULT,
complex*16, dimension( * )  WORK,
integer  NWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

ZDRVEV

Purpose:
    ZDRVEV  checks the nonsymmetric eigenvalue problem driver ZGEEV.

    When ZDRVEV 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, 7
    tests will be performed:

    (1)     | A * VR - VR * W | / ( n |A| ulp )

      Here VR is the matrix of unit right eigenvectors.
      W is a diagonal matrix with diagonal entries W(j).

    (2)     | A**H * VL - VL * W**H | / ( n |A| ulp )

      Here VL is the matrix of unit left eigenvectors, A**H is the
      conjugate-transpose of A, and W is as above.

    (3)     | |VR(i)| - 1 | / ulp and whether largest component real

      VR(i) denotes the i-th column of VR.

    (4)     | |VL(i)| - 1 | / ulp and whether largest component real

      VL(i) denotes the i-th column of VL.

    (5)     W(full) = W(partial)

      W(full) denotes the eigenvalues computed when both VR and VL
      are also computed, and W(partial) denotes the eigenvalues
      computed when only W, only W and VR, or only W and VL are
      computed.

    (6)     VR(full) = VR(partial)

      VR(full) denotes the right eigenvectors computed when both VR
      and VL are computed, and VR(partial) denotes the result
      when only VR is computed.

     (7)     VL(full) = VL(partial)

      VL(full) denotes the left eigenvectors computed when both VR
      and VL are also computed, and VL(partial) denotes the result
      when only VL is computed.

    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 unitary 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 |z| < 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,
          ZDRVEV 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, ZDRVEV
          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 ZDRVEV 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 ZGEEV.
[out]W
          W is COMPLEX*16 array, dimension (max(NN))
          The eigenvalues of A. On exit, W are the eigenvalues of
          the matrix in A.
[out]W1
          W1 is COMPLEX*16 array, dimension (max(NN))
          Like W, this array contains the eigenvalues of A,
          but those computed when ZGEEV only computes a partial
          eigendecomposition, i.e. not the eigenvalues and left
          and right eigenvectors.
[out]VL
          VL is COMPLEX*16 array, dimension (LDVL, max(NN))
          VL holds the computed left eigenvectors.
[in]LDVL
          LDVL is INTEGER
          Leading dimension of VL. Must be at least max(1,max(NN)).
[out]VR
          VR is COMPLEX*16 array, dimension (LDVR, max(NN))
          VR holds the computed right eigenvectors.
[in]LDVR
          LDVR is INTEGER
          Leading dimension of VR. Must be at least max(1,max(NN)).
[out]LRE
          LRE is COMPLEX*16 array, dimension (LDLRE, max(NN))
          LRE holds the computed right or left eigenvectors.
[in]LDLRE
          LDLRE is INTEGER
          Leading dimension of LRE. Must be at least max(1,max(NN)).
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (7)
          The values computed by the seven 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 (2*max(NN))
[out]IWORK
          IWORK is INTEGER 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) ).
          -14: LDVL < 1 or LDVL < NMAX, where NMAX is max( NN(j) ).
          -16: LDVR < 1 or LDVR < NMAX, where NMAX is max( NN(j) ).
          -18: LDLRE < 1 or LDLRE < NMAX, where NMAX is max( NN(j) ).
          -21: NWORK too small.
          If  ZLATMR, CLATMS, CLATME or ZGEEV 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)       Selectw 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 387 of file zdrvev.f.

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