LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sdrvev ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
real  THRESH,
integer  NOUNIT,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( lda, * )  H,
real, dimension( * )  WR,
real, dimension( * )  WI,
real, dimension( * )  WR1,
real, dimension( * )  WI1,
real, dimension( ldvl, * )  VL,
integer  LDVL,
real, dimension( ldvr, * )  VR,
integer  LDVR,
real, dimension( ldlre, * )  LRE,
integer  LDLRE,
real, dimension( 7 )  RESULT,
real, dimension( * )  WORK,
integer  NWORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

SDRVEV

Purpose:
    SDRVEV  checks the nonsymmetric eigenvalue problem driver SGEEV.

    When SDRVEV 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 block diagonal matrix, with a 1x1 block for each
      real eigenvalue and a 2x2 block for each complex conjugate
      pair.  If eigenvalues j and j+1 are a complex conjugate pair,
      so WR(j) = WR(j+1) = wr and WI(j) = - WI(j+1) = wi, then the
      2 x 2 block corresponding to the pair will be:

              (  wr  wi  )
              ( -wi  wr  )

      Such a block multiplying an n x 2 matrix  ( ur ui ) on the
      right will be the same as multiplying  ur + i*ui  by  wr + i*wi.

    (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 signs.
         (ULP = (first number larger than 1) - 1 )
    (5)  A diagonal matrix with geometrically spaced entries
         1, ..., ULP  and random signs.
    (6)  A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
         and random signs.

    (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 orthogonal and
         T has evenly spaced entries 1, ..., ULP with random signs
         on the diagonal and random O(1) entries in the upper
         triangle.

    (10) A matrix of the form  U' T U, where U is orthogonal and
         T has geometrically spaced entries 1, ..., ULP with random
         signs 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
         signs on the diagonal and random O(1) entries in the upper
         triangle.

    (12) A matrix of the form  U' T U, where U is orthogonal and
         T has real or complex conjugate paired eigenvalues randomly
         chosen from ( ULP, 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 signs 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 signs 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 signs 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 real or complex conjugate paired
         eigenvalues randomly chosen from ( ULP, 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,
          SDRVEV 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, SDRVEV
          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 SDRVEV 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 REAL 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 REAL array, dimension (LDA, max(NN))
          Another copy of the test matrix A, modified by SGEEV.
[out]WR
          WR is REAL array, dimension (max(NN))
[out]WI
          WI is REAL array, dimension (max(NN))

          The real and imaginary parts of the eigenvalues of A.
          On exit, WR + WI*i are the eigenvalues of the matrix in A.
[out]WR1
          WR1 is REAL array, dimension (max(NN))
[out]WI1
          WI1 is REAL array, dimension (max(NN))

          Like WR, WI, these arrays contain the eigenvalues of A,
          but those computed when SGEEV only computes a partial
          eigendecomposition, i.e. not the eigenvalues and left
          and right eigenvectors.
[out]VL
          VL is REAL 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 REAL 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 REAL 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 REAL 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]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) ).
          -16: LDVL < 1 or LDVL < NMAX, where NMAX is max( NN(j) ).
          -18: LDVR < 1 or LDVR < NMAX, where NMAX is max( NN(j) ).
          -20: LDLRE < 1 or LDLRE < NMAX, where NMAX is max( NN(j) ).
          -23: NWORK too small.
          If  SLATMR, SLATMS, SLATME or SGEEV 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.
Date
November 2011

Definition at line 408 of file sdrvev.f.

408 *
409 * -- LAPACK test routine (version 3.4.0) --
410 * -- LAPACK is a software package provided by Univ. of Tennessee, --
411 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
412 * November 2011
413 *
414 * .. Scalar Arguments ..
415  INTEGER info, lda, ldlre, ldvl, ldvr, nounit, nsizes,
416  $ ntypes, nwork
417  REAL thresh
418 * ..
419 * .. Array Arguments ..
420  LOGICAL dotype( * )
421  INTEGER iseed( 4 ), iwork( * ), nn( * )
422  REAL a( lda, * ), h( lda, * ), lre( ldlre, * ),
423  $ result( 7 ), vl( ldvl, * ), vr( ldvr, * ),
424  $ wi( * ), wi1( * ), work( * ), wr( * ), wr1( * )
425 * ..
426 *
427 * =====================================================================
428 *
429 * .. Parameters ..
430  REAL zero, one
431  parameter ( zero = 0.0e0, one = 1.0e0 )
432  REAL two
433  parameter ( two = 2.0e0 )
434  INTEGER maxtyp
435  parameter ( maxtyp = 21 )
436 * ..
437 * .. Local Scalars ..
438  LOGICAL badnn
439  CHARACTER*3 path
440  INTEGER iinfo, imode, itype, iwk, j, jcol, jj, jsize,
441  $ jtype, mtypes, n, nerrs, nfail, nmax,
442  $ nnwork, ntest, ntestf, ntestt
443  REAL anorm, cond, conds, ovfl, rtulp, rtulpi, tnrm,
444  $ ulp, ulpinv, unfl, vmx, vrmx, vtst
445 * ..
446 * .. Local Arrays ..
447  CHARACTER adumma( 1 )
448  INTEGER idumma( 1 ), ioldsd( 4 ), kconds( maxtyp ),
449  $ kmagn( maxtyp ), kmode( maxtyp ),
450  $ ktype( maxtyp )
451  REAL dum( 1 ), res( 2 )
452 * ..
453 * .. External Functions ..
454  REAL slamch, slapy2, snrm2
455  EXTERNAL slamch, slapy2, snrm2
456 * ..
457 * .. External Subroutines ..
458  EXTERNAL sgeev, sget22, slabad, slacpy, slasum, slatme,
460 * ..
461 * .. Intrinsic Functions ..
462  INTRINSIC abs, max, min, sqrt
463 * ..
464 * .. Data statements ..
465  DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
466  DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
467  $ 3, 1, 2, 3 /
468  DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
469  $ 1, 5, 5, 5, 4, 3, 1 /
470  DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
471 * ..
472 * .. Executable Statements ..
473 *
474  path( 1: 1 ) = 'Single precision'
475  path( 2: 3 ) = 'EV'
476 *
477 * Check for errors
478 *
479  ntestt = 0
480  ntestf = 0
481  info = 0
482 *
483 * Important constants
484 *
485  badnn = .false.
486  nmax = 0
487  DO 10 j = 1, nsizes
488  nmax = max( nmax, nn( j ) )
489  IF( nn( j ).LT.0 )
490  $ badnn = .true.
491  10 CONTINUE
492 *
493 * Check for errors
494 *
495  IF( nsizes.LT.0 ) THEN
496  info = -1
497  ELSE IF( badnn ) THEN
498  info = -2
499  ELSE IF( ntypes.LT.0 ) THEN
500  info = -3
501  ELSE IF( thresh.LT.zero ) THEN
502  info = -6
503  ELSE IF( nounit.LE.0 ) THEN
504  info = -7
505  ELSE IF( lda.LT.1 .OR. lda.LT.nmax ) THEN
506  info = -9
507  ELSE IF( ldvl.LT.1 .OR. ldvl.LT.nmax ) THEN
508  info = -16
509  ELSE IF( ldvr.LT.1 .OR. ldvr.LT.nmax ) THEN
510  info = -18
511  ELSE IF( ldlre.LT.1 .OR. ldlre.LT.nmax ) THEN
512  info = -20
513  ELSE IF( 5*nmax+2*nmax**2.GT.nwork ) THEN
514  info = -23
515  END IF
516 *
517  IF( info.NE.0 ) THEN
518  CALL xerbla( 'SDRVEV', -info )
519  RETURN
520  END IF
521 *
522 * Quick return if nothing to do
523 *
524  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
525  $ RETURN
526 *
527 * More Important constants
528 *
529  unfl = slamch( 'Safe minimum' )
530  ovfl = one / unfl
531  CALL slabad( unfl, ovfl )
532  ulp = slamch( 'Precision' )
533  ulpinv = one / ulp
534  rtulp = sqrt( ulp )
535  rtulpi = one / rtulp
536 *
537 * Loop over sizes, types
538 *
539  nerrs = 0
540 *
541  DO 270 jsize = 1, nsizes
542  n = nn( jsize )
543  IF( nsizes.NE.1 ) THEN
544  mtypes = min( maxtyp, ntypes )
545  ELSE
546  mtypes = min( maxtyp+1, ntypes )
547  END IF
548 *
549  DO 260 jtype = 1, mtypes
550  IF( .NOT.dotype( jtype ) )
551  $ GO TO 260
552 *
553 * Save ISEED in case of an error.
554 *
555  DO 20 j = 1, 4
556  ioldsd( j ) = iseed( j )
557  20 CONTINUE
558 *
559 * Compute "A"
560 *
561 * Control parameters:
562 *
563 * KMAGN KCONDS KMODE KTYPE
564 * =1 O(1) 1 clustered 1 zero
565 * =2 large large clustered 2 identity
566 * =3 small exponential Jordan
567 * =4 arithmetic diagonal, (w/ eigenvalues)
568 * =5 random log symmetric, w/ eigenvalues
569 * =6 random general, w/ eigenvalues
570 * =7 random diagonal
571 * =8 random symmetric
572 * =9 random general
573 * =10 random triangular
574 *
575  IF( mtypes.GT.maxtyp )
576  $ GO TO 90
577 *
578  itype = ktype( jtype )
579  imode = kmode( jtype )
580 *
581 * Compute norm
582 *
583  GO TO ( 30, 40, 50 )kmagn( jtype )
584 *
585  30 CONTINUE
586  anorm = one
587  GO TO 60
588 *
589  40 CONTINUE
590  anorm = ovfl*ulp
591  GO TO 60
592 *
593  50 CONTINUE
594  anorm = unfl*ulpinv
595  GO TO 60
596 *
597  60 CONTINUE
598 *
599  CALL slaset( 'Full', lda, n, zero, zero, a, lda )
600  iinfo = 0
601  cond = ulpinv
602 *
603 * Special Matrices -- Identity & Jordan block
604 *
605 * Zero
606 *
607  IF( itype.EQ.1 ) THEN
608  iinfo = 0
609 *
610  ELSE IF( itype.EQ.2 ) THEN
611 *
612 * Identity
613 *
614  DO 70 jcol = 1, n
615  a( jcol, jcol ) = anorm
616  70 CONTINUE
617 *
618  ELSE IF( itype.EQ.3 ) THEN
619 *
620 * Jordan Block
621 *
622  DO 80 jcol = 1, n
623  a( jcol, jcol ) = anorm
624  IF( jcol.GT.1 )
625  $ a( jcol, jcol-1 ) = one
626  80 CONTINUE
627 *
628  ELSE IF( itype.EQ.4 ) THEN
629 *
630 * Diagonal Matrix, [Eigen]values Specified
631 *
632  CALL slatms( n, n, 'S', iseed, 'S', work, imode, cond,
633  $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
634  $ iinfo )
635 *
636  ELSE IF( itype.EQ.5 ) THEN
637 *
638 * Symmetric, eigenvalues specified
639 *
640  CALL slatms( n, n, 'S', iseed, 'S', work, imode, cond,
641  $ anorm, n, n, 'N', a, lda, work( n+1 ),
642  $ iinfo )
643 *
644  ELSE IF( itype.EQ.6 ) THEN
645 *
646 * General, eigenvalues specified
647 *
648  IF( kconds( jtype ).EQ.1 ) THEN
649  conds = one
650  ELSE IF( kconds( jtype ).EQ.2 ) THEN
651  conds = rtulpi
652  ELSE
653  conds = zero
654  END IF
655 *
656  adumma( 1 ) = ' '
657  CALL slatme( n, 'S', iseed, work, imode, cond, one,
658  $ adumma, 'T', 'T', 'T', work( n+1 ), 4,
659  $ conds, n, n, anorm, a, lda, work( 2*n+1 ),
660  $ iinfo )
661 *
662  ELSE IF( itype.EQ.7 ) THEN
663 *
664 * Diagonal, random eigenvalues
665 *
666  CALL slatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
667  $ 'T', 'N', work( n+1 ), 1, one,
668  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
669  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
670 *
671  ELSE IF( itype.EQ.8 ) THEN
672 *
673 * Symmetric, random eigenvalues
674 *
675  CALL slatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
676  $ 'T', 'N', work( n+1 ), 1, one,
677  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
678  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
679 *
680  ELSE IF( itype.EQ.9 ) THEN
681 *
682 * General, random eigenvalues
683 *
684  CALL slatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
685  $ 'T', 'N', work( n+1 ), 1, one,
686  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
687  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
688  IF( n.GE.4 ) THEN
689  CALL slaset( 'Full', 2, n, zero, zero, a, lda )
690  CALL slaset( 'Full', n-3, 1, zero, zero, a( 3, 1 ),
691  $ lda )
692  CALL slaset( 'Full', n-3, 2, zero, zero, a( 3, n-1 ),
693  $ lda )
694  CALL slaset( 'Full', 1, n, zero, zero, a( n, 1 ),
695  $ lda )
696  END IF
697 *
698  ELSE IF( itype.EQ.10 ) THEN
699 *
700 * Triangular, random eigenvalues
701 *
702  CALL slatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
703  $ 'T', 'N', work( n+1 ), 1, one,
704  $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
705  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
706 *
707  ELSE
708 *
709  iinfo = 1
710  END IF
711 *
712  IF( iinfo.NE.0 ) THEN
713  WRITE( nounit, fmt = 9993 )'Generator', iinfo, n, jtype,
714  $ ioldsd
715  info = abs( iinfo )
716  RETURN
717  END IF
718 *
719  90 CONTINUE
720 *
721 * Test for minimal and generous workspace
722 *
723  DO 250 iwk = 1, 2
724  IF( iwk.EQ.1 ) THEN
725  nnwork = 4*n
726  ELSE
727  nnwork = 5*n + 2*n**2
728  END IF
729  nnwork = max( nnwork, 1 )
730 *
731 * Initialize RESULT
732 *
733  DO 100 j = 1, 7
734  result( j ) = -one
735  100 CONTINUE
736 *
737 * Compute eigenvalues and eigenvectors, and test them
738 *
739  CALL slacpy( 'F', n, n, a, lda, h, lda )
740  CALL sgeev( 'V', 'V', n, h, lda, wr, wi, vl, ldvl, vr,
741  $ ldvr, work, nnwork, iinfo )
742  IF( iinfo.NE.0 ) THEN
743  result( 1 ) = ulpinv
744  WRITE( nounit, fmt = 9993 )'SGEEV1', iinfo, n, jtype,
745  $ ioldsd
746  info = abs( iinfo )
747  GO TO 220
748  END IF
749 *
750 * Do Test (1)
751 *
752  CALL sget22( 'N', 'N', 'N', n, a, lda, vr, ldvr, wr, wi,
753  $ work, res )
754  result( 1 ) = res( 1 )
755 *
756 * Do Test (2)
757 *
758  CALL sget22( 'T', 'N', 'T', n, a, lda, vl, ldvl, wr, wi,
759  $ work, res )
760  result( 2 ) = res( 1 )
761 *
762 * Do Test (3)
763 *
764  DO 120 j = 1, n
765  tnrm = one
766  IF( wi( j ).EQ.zero ) THEN
767  tnrm = snrm2( n, vr( 1, j ), 1 )
768  ELSE IF( wi( j ).GT.zero ) THEN
769  tnrm = slapy2( snrm2( n, vr( 1, j ), 1 ),
770  $ snrm2( n, vr( 1, j+1 ), 1 ) )
771  END IF
772  result( 3 ) = max( result( 3 ),
773  $ min( ulpinv, abs( tnrm-one ) / ulp ) )
774  IF( wi( j ).GT.zero ) THEN
775  vmx = zero
776  vrmx = zero
777  DO 110 jj = 1, n
778  vtst = slapy2( vr( jj, j ), vr( jj, j+1 ) )
779  IF( vtst.GT.vmx )
780  $ vmx = vtst
781  IF( vr( jj, j+1 ).EQ.zero .AND.
782  $ abs( vr( jj, j ) ).GT.vrmx )
783  $ vrmx = abs( vr( jj, j ) )
784  110 CONTINUE
785  IF( vrmx / vmx.LT.one-two*ulp )
786  $ result( 3 ) = ulpinv
787  END IF
788  120 CONTINUE
789 *
790 * Do Test (4)
791 *
792  DO 140 j = 1, n
793  tnrm = one
794  IF( wi( j ).EQ.zero ) THEN
795  tnrm = snrm2( n, vl( 1, j ), 1 )
796  ELSE IF( wi( j ).GT.zero ) THEN
797  tnrm = slapy2( snrm2( n, vl( 1, j ), 1 ),
798  $ snrm2( n, vl( 1, j+1 ), 1 ) )
799  END IF
800  result( 4 ) = max( result( 4 ),
801  $ min( ulpinv, abs( tnrm-one ) / ulp ) )
802  IF( wi( j ).GT.zero ) THEN
803  vmx = zero
804  vrmx = zero
805  DO 130 jj = 1, n
806  vtst = slapy2( vl( jj, j ), vl( jj, j+1 ) )
807  IF( vtst.GT.vmx )
808  $ vmx = vtst
809  IF( vl( jj, j+1 ).EQ.zero .AND.
810  $ abs( vl( jj, j ) ).GT.vrmx )
811  $ vrmx = abs( vl( jj, j ) )
812  130 CONTINUE
813  IF( vrmx / vmx.LT.one-two*ulp )
814  $ result( 4 ) = ulpinv
815  END IF
816  140 CONTINUE
817 *
818 * Compute eigenvalues only, and test them
819 *
820  CALL slacpy( 'F', n, n, a, lda, h, lda )
821  CALL sgeev( 'N', 'N', n, h, lda, wr1, wi1, dum, 1, dum,
822  $ 1, work, nnwork, iinfo )
823  IF( iinfo.NE.0 ) THEN
824  result( 1 ) = ulpinv
825  WRITE( nounit, fmt = 9993 )'SGEEV2', iinfo, n, jtype,
826  $ ioldsd
827  info = abs( iinfo )
828  GO TO 220
829  END IF
830 *
831 * Do Test (5)
832 *
833  DO 150 j = 1, n
834  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
835  $ result( 5 ) = ulpinv
836  150 CONTINUE
837 *
838 * Compute eigenvalues and right eigenvectors, and test them
839 *
840  CALL slacpy( 'F', n, n, a, lda, h, lda )
841  CALL sgeev( 'N', 'V', n, h, lda, wr1, wi1, dum, 1, lre,
842  $ ldlre, work, nnwork, iinfo )
843  IF( iinfo.NE.0 ) THEN
844  result( 1 ) = ulpinv
845  WRITE( nounit, fmt = 9993 )'SGEEV3', iinfo, n, jtype,
846  $ ioldsd
847  info = abs( iinfo )
848  GO TO 220
849  END IF
850 *
851 * Do Test (5) again
852 *
853  DO 160 j = 1, n
854  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
855  $ result( 5 ) = ulpinv
856  160 CONTINUE
857 *
858 * Do Test (6)
859 *
860  DO 180 j = 1, n
861  DO 170 jj = 1, n
862  IF( vr( j, jj ).NE.lre( j, jj ) )
863  $ result( 6 ) = ulpinv
864  170 CONTINUE
865  180 CONTINUE
866 *
867 * Compute eigenvalues and left eigenvectors, and test them
868 *
869  CALL slacpy( 'F', n, n, a, lda, h, lda )
870  CALL sgeev( 'V', 'N', n, h, lda, wr1, wi1, lre, ldlre,
871  $ dum, 1, work, nnwork, iinfo )
872  IF( iinfo.NE.0 ) THEN
873  result( 1 ) = ulpinv
874  WRITE( nounit, fmt = 9993 )'SGEEV4', iinfo, n, jtype,
875  $ ioldsd
876  info = abs( iinfo )
877  GO TO 220
878  END IF
879 *
880 * Do Test (5) again
881 *
882  DO 190 j = 1, n
883  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
884  $ result( 5 ) = ulpinv
885  190 CONTINUE
886 *
887 * Do Test (7)
888 *
889  DO 210 j = 1, n
890  DO 200 jj = 1, n
891  IF( vl( j, jj ).NE.lre( j, jj ) )
892  $ result( 7 ) = ulpinv
893  200 CONTINUE
894  210 CONTINUE
895 *
896 * End of Loop -- Check for RESULT(j) > THRESH
897 *
898  220 CONTINUE
899 *
900  ntest = 0
901  nfail = 0
902  DO 230 j = 1, 7
903  IF( result( j ).GE.zero )
904  $ ntest = ntest + 1
905  IF( result( j ).GE.thresh )
906  $ nfail = nfail + 1
907  230 CONTINUE
908 *
909  IF( nfail.GT.0 )
910  $ ntestf = ntestf + 1
911  IF( ntestf.EQ.1 ) THEN
912  WRITE( nounit, fmt = 9999 )path
913  WRITE( nounit, fmt = 9998 )
914  WRITE( nounit, fmt = 9997 )
915  WRITE( nounit, fmt = 9996 )
916  WRITE( nounit, fmt = 9995 )thresh
917  ntestf = 2
918  END IF
919 *
920  DO 240 j = 1, 7
921  IF( result( j ).GE.thresh ) THEN
922  WRITE( nounit, fmt = 9994 )n, iwk, ioldsd, jtype,
923  $ j, result( j )
924  END IF
925  240 CONTINUE
926 *
927  nerrs = nerrs + nfail
928  ntestt = ntestt + ntest
929 *
930  250 CONTINUE
931  260 CONTINUE
932  270 CONTINUE
933 *
934 * Summary
935 *
936  CALL slasum( path, nounit, nerrs, ntestt )
937 *
938  9999 FORMAT( / 1x, a3, ' -- Real Eigenvalue-Eigenvector Decomposition',
939  $ ' Driver', / ' Matrix types (see SDRVEV for details): ' )
940 *
941  9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
942  $ ' ', ' 5=Diagonal: geometr. spaced entries.',
943  $ / ' 2=Identity matrix. ', ' 6=Diagona',
944  $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
945  $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
946  $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
947  $ 'mall, evenly spaced.' )
948  9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
949  $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
950  $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
951  $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
952  $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
953  $ 'lex ', / ' 12=Well-cond., random complex ', 6x, ' ',
954  $ ' 17=Ill-cond., large rand. complx ', / ' 13=Ill-condi',
955  $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
956  $ ' complx ' )
957  9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
958  $ 'with small random entries.', / ' 20=Matrix with large ran',
959  $ 'dom entries. ', / )
960  9995 FORMAT( ' Tests performed with test threshold =', f8.2,
961  $ / / ' 1 = | A VR - VR W | / ( n |A| ulp ) ',
962  $ / ' 2 = | transpose(A) VL - VL W | / ( n |A| ulp ) ',
963  $ / ' 3 = | |VR(i)| - 1 | / ulp ',
964  $ / ' 4 = | |VL(i)| - 1 | / ulp ',
965  $ / ' 5 = 0 if W same no matter if VR or VL computed,',
966  $ ' 1/ulp otherwise', /
967  $ ' 6 = 0 if VR same no matter if VL computed,',
968  $ ' 1/ulp otherwise', /
969  $ ' 7 = 0 if VL same no matter if VR computed,',
970  $ ' 1/ulp otherwise', / )
971  9994 FORMAT( ' N=', i5, ', IWK=', i2, ', seed=', 4( i4, ',' ),
972  $ ' type ', i2, ', test(', i2, ')=', g10.3 )
973  9993 FORMAT( ' SDRVEV: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
974  $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
975 *
976  RETURN
977 *
978 * End of SDRVEV
979 *
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:76
subroutine slatmr(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)
SLATMR
Definition: slatmr.f:473
real function slapy2(X, Y)
SLAPY2 returns sqrt(x2+y2).
Definition: slapy2.f:65
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
real function snrm2(N, X, INCX)
SNRM2
Definition: snrm2.f:56
subroutine sget22(TRANSA, TRANSE, TRANSW, N, A, LDA, E, LDE, WR, WI, WORK, RESULT)
SGET22
Definition: sget22.f:169
subroutine sgeev(JOBVL, JOBVR, N, A, LDA, WR, WI, VL, LDVL, VR, LDVR, WORK, LWORK, INFO)
SGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices ...
Definition: sgeev.f:193
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine slasum(TYPE, IOUNIT, IE, NRUN)
SLASUM
Definition: slasum.f:42
subroutine slatme(N, DIST, ISEED, D, MODE, COND, DMAX, EI, RSIGN, UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM, A, LDA, WORK, INFO)
SLATME
Definition: slatme.f:334

Here is the call graph for this function:

Here is the caller graph for this function: