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

◆ cdrvev()

subroutine cdrvev ( 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( * )  w,
complex, dimension( * )  w1,
complex, dimension( ldvl, * )  vl,
integer  ldvl,
complex, dimension( ldvr, * )  vr,
integer  ldvr,
complex, dimension( ldlre, * )  lre,
integer  ldlre,
real, dimension( 7 )  result,
complex, dimension( * )  work,
integer  nwork,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  info 
)

CDRVEV

Purpose:
    CDRVEV  checks the nonsymmetric eigenvalue problem driver CGEEV.

    When CDRVEV 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,
          CDRVEV 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, CDRVEV
          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 CDRVEV 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 CGEEV.
[out]W
          W is COMPLEX array, dimension (max(NN))
          The eigenvalues of A. On exit, W are the eigenvalues of
          the matrix in A.
[out]W1
          W1 is COMPLEX array, dimension (max(NN))
          Like W, this array contains the eigenvalues of A,
          but those computed when CGEEV only computes a partial
          eigendecomposition, i.e. not the eigenvalues and left
          and right eigenvectors.
[out]VL
          VL is COMPLEX 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 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 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 REAL 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 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 (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  CLATMR, CLATMS, CLATME or CGEEV 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 cdrvev.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 REAL THRESH
400* ..
401* .. Array Arguments ..
402 LOGICAL DOTYPE( * )
403 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
404 REAL RESULT( 7 ), RWORK( * )
405 COMPLEX A( LDA, * ), H( LDA, * ), LRE( LDLRE, * ),
406 $ VL( LDVL, * ), VR( LDVR, * ), W( * ), W1( * ),
407 $ WORK( * )
408* ..
409*
410* =====================================================================
411*
412* .. Parameters ..
413 COMPLEX CZERO
414 parameter( czero = ( 0.0e+0, 0.0e+0 ) )
415 COMPLEX CONE
416 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
417 REAL ZERO, ONE
418 parameter( zero = 0.0e+0, one = 1.0e+0 )
419 REAL TWO
420 parameter( two = 2.0e+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,
429 $ NNWORK, NTEST, NTESTF, NTESTT
430 REAL 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 REAL RES( 2 )
438 COMPLEX DUM( 1 )
439* ..
440* .. External Functions ..
441 REAL SCNRM2, SLAMCH
442 EXTERNAL scnrm2, slamch
443* ..
444* .. External Subroutines ..
445 EXTERNAL cgeev, cget22, clacpy, clatme, clatmr, clatms,
447* ..
448* .. Intrinsic Functions ..
449 INTRINSIC abs, aimag, cmplx, max, min, real, 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 ) = 'Complex 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( 'CDRVEV', -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 = slamch( 'Safe minimum' )
517 ovfl = one / unfl
518 ulp = slamch( 'Precision' )
519 ulpinv = one / ulp
520 rtulp = sqrt( ulp )
521 rtulpi = one / rtulp
522*
523* Loop over sizes, types
524*
525 nerrs = 0
526*
527 DO 270 jsize = 1, nsizes
528 n = nn( jsize )
529 IF( nsizes.NE.1 ) THEN
530 mtypes = min( maxtyp, ntypes )
531 ELSE
532 mtypes = min( maxtyp+1, ntypes )
533 END IF
534*
535 DO 260 jtype = 1, mtypes
536 IF( .NOT.dotype( jtype ) )
537 $ GO TO 260
538*
539* Save ISEED in case of an error.
540*
541 DO 20 j = 1, 4
542 ioldsd( j ) = iseed( j )
543 20 CONTINUE
544*
545* Compute "A"
546*
547* Control parameters:
548*
549* KMAGN KCONDS KMODE KTYPE
550* =1 O(1) 1 clustered 1 zero
551* =2 large large clustered 2 identity
552* =3 small exponential Jordan
553* =4 arithmetic diagonal, (w/ eigenvalues)
554* =5 random log symmetric, w/ eigenvalues
555* =6 random general, w/ eigenvalues
556* =7 random diagonal
557* =8 random symmetric
558* =9 random general
559* =10 random triangular
560*
561 IF( mtypes.GT.maxtyp )
562 $ GO TO 90
563*
564 itype = ktype( jtype )
565 imode = kmode( jtype )
566*
567* Compute norm
568*
569 GO TO ( 30, 40, 50 )kmagn( jtype )
570*
571 30 CONTINUE
572 anorm = one
573 GO TO 60
574*
575 40 CONTINUE
576 anorm = ovfl*ulp
577 GO TO 60
578*
579 50 CONTINUE
580 anorm = unfl*ulpinv
581 GO TO 60
582*
583 60 CONTINUE
584*
585 CALL claset( 'Full', lda, n, czero, czero, a, lda )
586 iinfo = 0
587 cond = ulpinv
588*
589* Special Matrices -- Identity & Jordan block
590*
591* Zero
592*
593 IF( itype.EQ.1 ) THEN
594 iinfo = 0
595*
596 ELSE IF( itype.EQ.2 ) THEN
597*
598* Identity
599*
600 DO 70 jcol = 1, n
601 a( jcol, jcol ) = cmplx( anorm )
602 70 CONTINUE
603*
604 ELSE IF( itype.EQ.3 ) THEN
605*
606* Jordan Block
607*
608 DO 80 jcol = 1, n
609 a( jcol, jcol ) = cmplx( anorm )
610 IF( jcol.GT.1 )
611 $ a( jcol, jcol-1 ) = cone
612 80 CONTINUE
613*
614 ELSE IF( itype.EQ.4 ) THEN
615*
616* Diagonal Matrix, [Eigen]values Specified
617*
618 CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
619 $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
620 $ iinfo )
621*
622 ELSE IF( itype.EQ.5 ) THEN
623*
624* Hermitian, eigenvalues specified
625*
626 CALL clatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
627 $ anorm, n, n, 'N', a, lda, work( n+1 ),
628 $ iinfo )
629*
630 ELSE IF( itype.EQ.6 ) THEN
631*
632* General, eigenvalues specified
633*
634 IF( kconds( jtype ).EQ.1 ) THEN
635 conds = one
636 ELSE IF( kconds( jtype ).EQ.2 ) THEN
637 conds = rtulpi
638 ELSE
639 conds = zero
640 END IF
641*
642 CALL clatme( n, 'D', iseed, work, imode, cond, cone,
643 $ 'T', 'T', 'T', rwork, 4, conds, n, n,
644 $ anorm, a, lda, work( 2*n+1 ), iinfo )
645*
646 ELSE IF( itype.EQ.7 ) THEN
647*
648* Diagonal, random eigenvalues
649*
650 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
651 $ 'T', 'N', work( n+1 ), 1, one,
652 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
653 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
654*
655 ELSE IF( itype.EQ.8 ) THEN
656*
657* Symmetric, random eigenvalues
658*
659 CALL clatmr( n, n, 'D', iseed, 'H', work, 6, one, cone,
660 $ 'T', 'N', work( n+1 ), 1, one,
661 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
662 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
663*
664 ELSE IF( itype.EQ.9 ) THEN
665*
666* General, random eigenvalues
667*
668 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
669 $ 'T', 'N', work( n+1 ), 1, one,
670 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
671 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
672 IF( n.GE.4 ) THEN
673 CALL claset( 'Full', 2, n, czero, czero, a, lda )
674 CALL claset( 'Full', n-3, 1, czero, czero, a( 3, 1 ),
675 $ lda )
676 CALL claset( 'Full', n-3, 2, czero, czero,
677 $ a( 3, n-1 ), lda )
678 CALL claset( 'Full', 1, n, czero, czero, a( n, 1 ),
679 $ lda )
680 END IF
681*
682 ELSE IF( itype.EQ.10 ) THEN
683*
684* Triangular, random eigenvalues
685*
686 CALL clatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
687 $ 'T', 'N', work( n+1 ), 1, one,
688 $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
689 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
690*
691 ELSE
692*
693 iinfo = 1
694 END IF
695*
696 IF( iinfo.NE.0 ) THEN
697 WRITE( nounit, fmt = 9993 )'Generator', iinfo, n, jtype,
698 $ ioldsd
699 info = abs( iinfo )
700 RETURN
701 END IF
702*
703 90 CONTINUE
704*
705* Test for minimal and generous workspace
706*
707 DO 250 iwk = 1, 2
708 IF( iwk.EQ.1 ) THEN
709 nnwork = 2*n
710 ELSE
711 nnwork = 5*n + 2*n**2
712 END IF
713 nnwork = max( nnwork, 1 )
714*
715* Initialize RESULT
716*
717 DO 100 j = 1, 7
718 result( j ) = -one
719 100 CONTINUE
720*
721* Compute eigenvalues and eigenvectors, and test them
722*
723 CALL clacpy( 'F', n, n, a, lda, h, lda )
724 CALL cgeev( 'V', 'V', n, h, lda, w, vl, ldvl, vr, ldvr,
725 $ work, nnwork, rwork, iinfo )
726 IF( iinfo.NE.0 ) THEN
727 result( 1 ) = ulpinv
728 WRITE( nounit, fmt = 9993 )'CGEEV1', iinfo, n, jtype,
729 $ ioldsd
730 info = abs( iinfo )
731 GO TO 220
732 END IF
733*
734* Do Test (1)
735*
736 CALL cget22( 'N', 'N', 'N', n, a, lda, vr, ldvr, w, work,
737 $ rwork, res )
738 result( 1 ) = res( 1 )
739*
740* Do Test (2)
741*
742 CALL cget22( 'C', 'N', 'C', n, a, lda, vl, ldvl, w, work,
743 $ rwork, res )
744 result( 2 ) = res( 1 )
745*
746* Do Test (3)
747*
748 DO 120 j = 1, n
749 tnrm = scnrm2( n, vr( 1, j ), 1 )
750 result( 3 ) = max( result( 3 ),
751 $ min( ulpinv, abs( tnrm-one ) / ulp ) )
752 vmx = zero
753 vrmx = zero
754 DO 110 jj = 1, n
755 vtst = abs( vr( jj, j ) )
756 IF( vtst.GT.vmx )
757 $ vmx = vtst
758 IF( aimag( vr( jj, j ) ).EQ.zero .AND.
759 $ abs( real( vr( jj, j ) ) ).GT.vrmx )
760 $ vrmx = abs( real( vr( jj, j ) ) )
761 110 CONTINUE
762 IF( vrmx / vmx.LT.one-two*ulp )
763 $ result( 3 ) = ulpinv
764 120 CONTINUE
765*
766* Do Test (4)
767*
768 DO 140 j = 1, n
769 tnrm = scnrm2( n, vl( 1, j ), 1 )
770 result( 4 ) = max( result( 4 ),
771 $ min( ulpinv, abs( tnrm-one ) / ulp ) )
772 vmx = zero
773 vrmx = zero
774 DO 130 jj = 1, n
775 vtst = abs( vl( jj, j ) )
776 IF( vtst.GT.vmx )
777 $ vmx = vtst
778 IF( aimag( vl( jj, j ) ).EQ.zero .AND.
779 $ abs( real( vl( jj, j ) ) ).GT.vrmx )
780 $ vrmx = abs( real( vl( jj, j ) ) )
781 130 CONTINUE
782 IF( vrmx / vmx.LT.one-two*ulp )
783 $ result( 4 ) = ulpinv
784 140 CONTINUE
785*
786* Compute eigenvalues only, and test them
787*
788 CALL clacpy( 'F', n, n, a, lda, h, lda )
789 CALL cgeev( 'N', 'N', n, h, lda, w1, dum, 1, dum, 1,
790 $ work, nnwork, rwork, iinfo )
791 IF( iinfo.NE.0 ) THEN
792 result( 1 ) = ulpinv
793 WRITE( nounit, fmt = 9993 )'CGEEV2', iinfo, n, jtype,
794 $ ioldsd
795 info = abs( iinfo )
796 GO TO 220
797 END IF
798*
799* Do Test (5)
800*
801 DO 150 j = 1, n
802 IF( w( j ).NE.w1( j ) )
803 $ result( 5 ) = ulpinv
804 150 CONTINUE
805*
806* Compute eigenvalues and right eigenvectors, and test them
807*
808 CALL clacpy( 'F', n, n, a, lda, h, lda )
809 CALL cgeev( 'N', 'V', n, h, lda, w1, dum, 1, lre, ldlre,
810 $ work, nnwork, rwork, iinfo )
811 IF( iinfo.NE.0 ) THEN
812 result( 1 ) = ulpinv
813 WRITE( nounit, fmt = 9993 )'CGEEV3', iinfo, n, jtype,
814 $ ioldsd
815 info = abs( iinfo )
816 GO TO 220
817 END IF
818*
819* Do Test (5) again
820*
821 DO 160 j = 1, n
822 IF( w( j ).NE.w1( j ) )
823 $ result( 5 ) = ulpinv
824 160 CONTINUE
825*
826* Do Test (6)
827*
828 DO 180 j = 1, n
829 DO 170 jj = 1, n
830 IF( vr( j, jj ).NE.lre( j, jj ) )
831 $ result( 6 ) = ulpinv
832 170 CONTINUE
833 180 CONTINUE
834*
835* Compute eigenvalues and left eigenvectors, and test them
836*
837 CALL clacpy( 'F', n, n, a, lda, h, lda )
838 CALL cgeev( 'V', 'N', n, h, lda, w1, lre, ldlre, dum, 1,
839 $ work, nnwork, rwork, iinfo )
840 IF( iinfo.NE.0 ) THEN
841 result( 1 ) = ulpinv
842 WRITE( nounit, fmt = 9993 )'CGEEV4', iinfo, n, jtype,
843 $ ioldsd
844 info = abs( iinfo )
845 GO TO 220
846 END IF
847*
848* Do Test (5) again
849*
850 DO 190 j = 1, n
851 IF( w( j ).NE.w1( j ) )
852 $ result( 5 ) = ulpinv
853 190 CONTINUE
854*
855* Do Test (7)
856*
857 DO 210 j = 1, n
858 DO 200 jj = 1, n
859 IF( vl( j, jj ).NE.lre( j, jj ) )
860 $ result( 7 ) = ulpinv
861 200 CONTINUE
862 210 CONTINUE
863*
864* End of Loop -- Check for RESULT(j) > THRESH
865*
866 220 CONTINUE
867*
868 ntest = 0
869 nfail = 0
870 DO 230 j = 1, 7
871 IF( result( j ).GE.zero )
872 $ ntest = ntest + 1
873 IF( result( j ).GE.thresh )
874 $ nfail = nfail + 1
875 230 CONTINUE
876*
877 IF( nfail.GT.0 )
878 $ ntestf = ntestf + 1
879 IF( ntestf.EQ.1 ) THEN
880 WRITE( nounit, fmt = 9999 )path
881 WRITE( nounit, fmt = 9998 )
882 WRITE( nounit, fmt = 9997 )
883 WRITE( nounit, fmt = 9996 )
884 WRITE( nounit, fmt = 9995 )thresh
885 ntestf = 2
886 END IF
887*
888 DO 240 j = 1, 7
889 IF( result( j ).GE.thresh ) THEN
890 WRITE( nounit, fmt = 9994 )n, iwk, ioldsd, jtype,
891 $ j, result( j )
892 END IF
893 240 CONTINUE
894*
895 nerrs = nerrs + nfail
896 ntestt = ntestt + ntest
897*
898 250 CONTINUE
899 260 CONTINUE
900 270 CONTINUE
901*
902* Summary
903*
904 CALL slasum( path, nounit, nerrs, ntestt )
905*
906 9999 FORMAT( / 1x, a3, ' -- Complex Eigenvalue-Eigenvector ',
907 $ 'Decomposition Driver', /
908 $ ' Matrix types (see CDRVEV for details): ' )
909*
910 9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
911 $ ' ', ' 5=Diagonal: geometr. spaced entries.',
912 $ / ' 2=Identity matrix. ', ' 6=Diagona',
913 $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
914 $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
915 $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
916 $ 'mall, evenly spaced.' )
917 9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
918 $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
919 $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
920 $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
921 $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
922 $ 'lex ', a6, / ' 12=Well-cond., random complex ', a6, ' ',
923 $ ' 17=Ill-cond., large rand. complx ', a4, / ' 13=Ill-condi',
924 $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
925 $ ' complx ', a4 )
926 9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
927 $ 'with small random entries.', / ' 20=Matrix with large ran',
928 $ 'dom entries. ', / )
929 9995 FORMAT( ' Tests performed with test threshold =', f8.2,
930 $ / / ' 1 = | A VR - VR W | / ( n |A| ulp ) ',
931 $ / ' 2 = | conj-trans(A) VL - VL conj-trans(W) | /',
932 $ ' ( n |A| ulp ) ', / ' 3 = | |VR(i)| - 1 | / ulp ',
933 $ / ' 4 = | |VL(i)| - 1 | / ulp ',
934 $ / ' 5 = 0 if W same no matter if VR or VL computed,',
935 $ ' 1/ulp otherwise', /
936 $ ' 6 = 0 if VR same no matter if VL computed,',
937 $ ' 1/ulp otherwise', /
938 $ ' 7 = 0 if VL same no matter if VR computed,',
939 $ ' 1/ulp otherwise', / )
940 9994 FORMAT( ' N=', i5, ', IWK=', i2, ', seed=', 4( i4, ',' ),
941 $ ' type ', i2, ', test(', i2, ')=', g10.3 )
942 9993 FORMAT( ' CDRVEV: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
943 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
944*
945 RETURN
946*
947* End of CDRVEV
948*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cget22(transa, transe, transw, n, a, lda, e, lde, w, work, rwork, result)
CGET22
Definition cget22.f:144
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
subroutine cgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
CGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
Definition cgeev.f:180
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
real(wp) function scnrm2(n, x, incx)
SCNRM2
Definition scnrm2.f90:90
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: