LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrges ( 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, * )  B,
complex*16, dimension( lda, * )  S,
complex*16, dimension( lda, * )  T,
complex*16, dimension( ldq, * )  Q,
integer  LDQ,
complex*16, dimension( ldq, * )  Z,
complex*16, dimension( * )  ALPHA,
complex*16, dimension( * )  BETA,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
double precision, dimension( 13 )  RESULT,
logical, dimension( * )  BWORK,
integer  INFO 
)

ZDRGES

Purpose:
 ZDRGES checks the nonsymmetric generalized eigenvalue (Schur form)
 problem driver ZGGES.

 ZGGES factors A and B as Q*S*Z'  and Q*T*Z' , where ' means conjugate
 transpose, S and T are  upper triangular (i.e., in generalized Schur
 form), and Q and Z are unitary. It also computes the generalized
 eigenvalues (alpha(j),beta(j)), j=1,...,n.  Thus,
 w(j) = alpha(j)/beta(j) is a root of the characteristic equation

                 det( A - w(j) B ) = 0

 Optionally it also reorder the eigenvalues so that a selected
 cluster of eigenvalues appears in the leading diagonal block of the
 Schur forms.

 When ZDRGES 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, a pair of matrices (A, B) will be generated
 and used for testing. For each matrix pair, the following 13 tests
 will be performed and compared with the threshold THRESH except
 the tests (5), (11) and (13).


 (1)   | A - Q S Z' | / ( |A| n ulp ) (no sorting of eigenvalues)


 (2)   | B - Q T Z' | / ( |B| n ulp ) (no sorting of eigenvalues)


 (3)   | I - QQ' | / ( n ulp ) (no sorting of eigenvalues)


 (4)   | I - ZZ' | / ( n ulp ) (no sorting of eigenvalues)

 (5)   if A is in Schur form (i.e. triangular form) (no sorting of
       eigenvalues)

 (6)   if eigenvalues = diagonal elements of the Schur form (S, T),
       i.e., test the maximum over j of D(j)  where:

                     |alpha(j) - S(j,j)|        |beta(j) - T(j,j)|
           D(j) = ------------------------ + -----------------------
                  max(|alpha(j)|,|S(j,j)|)   max(|beta(j)|,|T(j,j)|)

       (no sorting of eigenvalues)

 (7)   | (A,B) - Q (S,T) Z' | / ( |(A,B)| n ulp )
       (with sorting of eigenvalues).

 (8)   | I - QQ' | / ( n ulp ) (with sorting of eigenvalues).

 (9)   | I - ZZ' | / ( n ulp ) (with sorting of eigenvalues).

 (10)  if A is in Schur form (i.e. quasi-triangular form)
       (with sorting of eigenvalues).

 (11)  if eigenvalues = diagonal elements of the Schur form (S, T),
       i.e. test the maximum over j of D(j)  where:

                     |alpha(j) - S(j,j)|        |beta(j) - T(j,j)|
           D(j) = ------------------------ + -----------------------
                  max(|alpha(j)|,|S(j,j)|)   max(|beta(j)|,|T(j,j)|)

       (with sorting of eigenvalues).

 (12)  if sorting worked and SDIM is the number of eigenvalues
       which were CELECTed.

 Test Matrices
 =============

 The sizes of the test matrices 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)  ( 0, 0 )         (a pair of zero matrices)

 (2)  ( I, 0 )         (an identity and a zero matrix)

 (3)  ( 0, I )         (an identity and a zero matrix)

 (4)  ( I, I )         (a pair of identity matrices)

         t   t
 (5)  ( J , J  )       (a pair of transposed Jordan blocks)

                                     t                ( I   0  )
 (6)  ( X, Y )         where  X = ( J   0  )  and Y = (      t )
                                  ( 0   I  )          ( 0   J  )
                       and I is a k x k identity and J a (k+1)x(k+1)
                       Jordan block; k=(N-1)/2

 (7)  ( D, I )         where D is diag( 0, 1,..., N-1 ) (a diagonal
                       matrix with those diagonal entries.)
 (8)  ( I, D )

 (9)  ( big*D, small*I ) where "big" is near overflow and small=1/big

 (10) ( small*D, big*I )

 (11) ( big*I, small*D )

 (12) ( small*I, big*D )

 (13) ( big*D, big*I )

 (14) ( small*D, small*I )

 (15) ( D1, D2 )        where D1 is diag( 0, 0, 1, ..., N-3, 0 ) and
                        D2 is diag( 0, N-3, N-4,..., 1, 0, 0 )
           t   t
 (16) Q ( J , J ) Z     where Q and Z are random orthogonal matrices.

 (17) Q ( T1, T2 ) Z    where T1 and T2 are upper triangular matrices
                        with random O(1) entries above the diagonal
                        and diagonal entries diag(T1) =
                        ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
                        ( 0, N-3, N-4,..., 1, 0, 0 )

 (18) Q ( T1, T2 ) Z    diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
                        diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
                        s = machine precision.

 (19) Q ( T1, T2 ) Z    diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 )
                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 )

                                                        N-5
 (20) Q ( T1, T2 ) Z    diag(T1)=( 0, 0, 1, 1, a, ..., a   =s, 0 )
                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )

 (21) Q ( T1, T2 ) Z    diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 )
                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
                        where r1,..., r(N-4) are random.

 (22) Q ( big*T1, small*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (23) Q ( small*T1, big*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (24) Q ( small*T1, small*T2 ) Z  diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (25) Q ( big*T1, big*T2 ) Z      diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (26) Q ( T1, T2 ) Z     where T1 and T2 are random upper-triangular
                         matrices.
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of sizes of matrices to use.  If it is zero,
          DDRGES does nothing.  NSIZES >= 0.
[in]NN
          NN is INTEGER array, dimension (NSIZES)
          An array containing the sizes to be used for the matrices.
          Zero values will be skipped.  NN >= 0.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, DDRGES
          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 on input.
          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 DDRGES 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.  THRESH >= 0.
[in]NOUNIT
          NOUNIT is INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
[in,out]A
          A is COMPLEX*16 array, dimension(LDA, max(NN))
          Used to hold the original A matrix.  Used as input only
          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
          DOTYPE(MAXTYP+1)=.TRUE.
[in]LDA
          LDA is INTEGER
          The leading dimension of A, B, S, and T.
          It must be at least 1 and at least max( NN ).
[in,out]B
          B is COMPLEX*16 array, dimension(LDA, max(NN))
          Used to hold the original B matrix.  Used as input only
          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
          DOTYPE(MAXTYP+1)=.TRUE.
[out]S
          S is COMPLEX*16 array, dimension (LDA, max(NN))
          The Schur form matrix computed from A by ZGGES.  On exit, S
          contains the Schur form matrix corresponding to the matrix
          in A.
[out]T
          T is COMPLEX*16 array, dimension (LDA, max(NN))
          The upper triangular matrix computed from B by ZGGES.
[out]Q
          Q is COMPLEX*16 array, dimension (LDQ, max(NN))
          The (left) orthogonal matrix computed by ZGGES.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of Q and Z. It must
          be at least 1 and at least max( NN ).
[out]Z
          Z is COMPLEX*16 array, dimension( LDQ, max(NN) )
          The (right) orthogonal matrix computed by ZGGES.
[out]ALPHA
          ALPHA is COMPLEX*16 array, dimension (max(NN))
[out]BETA
          BETA is COMPLEX*16 array, dimension (max(NN))

          The generalized eigenvalues of (A,B) computed by ZGGES.
          ALPHA(k) / BETA(k) is the k-th generalized eigenvalue of A
          and B.
[out]WORK
          WORK is COMPLEX*16 array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= 3*N*N.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension ( 8*N )
          Real workspace.
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (15)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid overflow.
[out]BWORK
          BWORK is LOGICAL array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  A routine returned an error code.  INFO is the
                absolute value of the INFO value returned.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
June 2016

Definition at line 383 of file zdrges.f.

383 *
384 * -- LAPACK test routine (version 3.6.1) --
385 * -- LAPACK is a software package provided by Univ. of Tennessee, --
386 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
387 * June 2016
388 *
389 * .. Scalar Arguments ..
390  INTEGER info, lda, ldq, lwork, nounit, nsizes, ntypes
391  DOUBLE PRECISION thresh
392 * ..
393 * .. Array Arguments ..
394  LOGICAL bwork( * ), dotype( * )
395  INTEGER iseed( 4 ), nn( * )
396  DOUBLE PRECISION result( 13 ), rwork( * )
397  COMPLEX*16 a( lda, * ), alpha( * ), b( lda, * ),
398  $ beta( * ), q( ldq, * ), s( lda, * ),
399  $ t( lda, * ), work( * ), z( ldq, * )
400 * ..
401 *
402 * =====================================================================
403 *
404 * .. Parameters ..
405  DOUBLE PRECISION zero, one
406  parameter ( zero = 0.0d+0, one = 1.0d+0 )
407  COMPLEX*16 czero, cone
408  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
409  $ cone = ( 1.0d+0, 0.0d+0 ) )
410  INTEGER maxtyp
411  parameter ( maxtyp = 26 )
412 * ..
413 * .. Local Scalars ..
414  LOGICAL badnn, ilabad
415  CHARACTER sort
416  INTEGER i, iadd, iinfo, in, isort, j, jc, jr, jsize,
417  $ jtype, knteig, maxwrk, minwrk, mtypes, n, n1,
418  $ nb, nerrs, nmats, nmax, ntest, ntestt, rsub,
419  $ sdim
420  DOUBLE PRECISION safmax, safmin, temp1, temp2, ulp, ulpinv
421  COMPLEX*16 ctemp, x
422 * ..
423 * .. Local Arrays ..
424  LOGICAL lasign( maxtyp ), lbsign( maxtyp )
425  INTEGER ioldsd( 4 ), kadd( 6 ), kamagn( maxtyp ),
426  $ katype( maxtyp ), kazero( maxtyp ),
427  $ kbmagn( maxtyp ), kbtype( maxtyp ),
428  $ kbzero( maxtyp ), kclass( maxtyp ),
429  $ ktrian( maxtyp ), kz1( 6 ), kz2( 6 )
430  DOUBLE PRECISION rmagn( 0: 3 )
431 * ..
432 * .. External Functions ..
433  LOGICAL zlctes
434  INTEGER ilaenv
435  DOUBLE PRECISION dlamch
436  COMPLEX*16 zlarnd
437  EXTERNAL zlctes, ilaenv, dlamch, zlarnd
438 * ..
439 * .. External Subroutines ..
440  EXTERNAL alasvm, dlabad, xerbla, zget51, zget54, zgges,
442 * ..
443 * .. Intrinsic Functions ..
444  INTRINSIC abs, dble, dconjg, dimag, max, min, sign
445 * ..
446 * .. Statement Functions ..
447  DOUBLE PRECISION abs1
448 * ..
449 * .. Statement Function definitions ..
450  abs1( x ) = abs( dble( x ) ) + abs( dimag( x ) )
451 * ..
452 * .. Data statements ..
453  DATA kclass / 15*1, 10*2, 1*3 /
454  DATA kz1 / 0, 1, 2, 1, 3, 3 /
455  DATA kz2 / 0, 0, 1, 2, 1, 1 /
456  DATA kadd / 0, 0, 0, 0, 3, 2 /
457  DATA katype / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
458  $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 /
459  DATA kbtype / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
460  $ 1, 1, -4, 2, -4, 8*8, 0 /
461  DATA kazero / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
462  $ 4*5, 4*3, 1 /
463  DATA kbzero / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
464  $ 4*6, 4*4, 1 /
465  DATA kamagn / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
466  $ 2, 1 /
467  DATA kbmagn / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
468  $ 2, 1 /
469  DATA ktrian / 16*0, 10*1 /
470  DATA lasign / 6*.false., .true., .false., 2*.true.,
471  $ 2*.false., 3*.true., .false., .true.,
472  $ 3*.false., 5*.true., .false. /
473  DATA lbsign / 7*.false., .true., 2*.false.,
474  $ 2*.true., 2*.false., .true., .false., .true.,
475  $ 9*.false. /
476 * ..
477 * .. Executable Statements ..
478 *
479 * Check for errors
480 *
481  info = 0
482 *
483  badnn = .false.
484  nmax = 1
485  DO 10 j = 1, nsizes
486  nmax = max( nmax, nn( j ) )
487  IF( nn( j ).LT.0 )
488  $ badnn = .true.
489  10 CONTINUE
490 *
491  IF( nsizes.LT.0 ) THEN
492  info = -1
493  ELSE IF( badnn ) THEN
494  info = -2
495  ELSE IF( ntypes.LT.0 ) THEN
496  info = -3
497  ELSE IF( thresh.LT.zero ) THEN
498  info = -6
499  ELSE IF( lda.LE.1 .OR. lda.LT.nmax ) THEN
500  info = -9
501  ELSE IF( ldq.LE.1 .OR. ldq.LT.nmax ) THEN
502  info = -14
503  END IF
504 *
505 * Compute workspace
506 * (Note: Comments in the code beginning "Workspace:" describe the
507 * minimal amount of workspace needed at that point in the code,
508 * as well as the preferred amount for good performance.
509 * NB refers to the optimal block size for the immediately
510 * following subroutine, as returned by ILAENV.
511 *
512  minwrk = 1
513  IF( info.EQ.0 .AND. lwork.GE.1 ) THEN
514  minwrk = 3*nmax*nmax
515  nb = max( 1, ilaenv( 1, 'ZGEQRF', ' ', nmax, nmax, -1, -1 ),
516  $ ilaenv( 1, 'ZUNMQR', 'LC', nmax, nmax, nmax, -1 ),
517  $ ilaenv( 1, 'ZUNGQR', ' ', nmax, nmax, nmax, -1 ) )
518  maxwrk = max( nmax+nmax*nb, 3*nmax*nmax )
519  work( 1 ) = maxwrk
520  END IF
521 *
522  IF( lwork.LT.minwrk )
523  $ info = -19
524 *
525  IF( info.NE.0 ) THEN
526  CALL xerbla( 'ZDRGES', -info )
527  RETURN
528  END IF
529 *
530 * Quick return if possible
531 *
532  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
533  $ RETURN
534 *
535  ulp = dlamch( 'Precision' )
536  safmin = dlamch( 'Safe minimum' )
537  safmin = safmin / ulp
538  safmax = one / safmin
539  CALL dlabad( safmin, safmax )
540  ulpinv = one / ulp
541 *
542 * The values RMAGN(2:3) depend on N, see below.
543 *
544  rmagn( 0 ) = zero
545  rmagn( 1 ) = one
546 *
547 * Loop over matrix sizes
548 *
549  ntestt = 0
550  nerrs = 0
551  nmats = 0
552 *
553  DO 190 jsize = 1, nsizes
554  n = nn( jsize )
555  n1 = max( 1, n )
556  rmagn( 2 ) = safmax*ulp / dble( n1 )
557  rmagn( 3 ) = safmin*ulpinv*dble( n1 )
558 *
559  IF( nsizes.NE.1 ) THEN
560  mtypes = min( maxtyp, ntypes )
561  ELSE
562  mtypes = min( maxtyp+1, ntypes )
563  END IF
564 *
565 * Loop over matrix types
566 *
567  DO 180 jtype = 1, mtypes
568  IF( .NOT.dotype( jtype ) )
569  $ GO TO 180
570  nmats = nmats + 1
571  ntest = 0
572 *
573 * Save ISEED in case of an error.
574 *
575  DO 20 j = 1, 4
576  ioldsd( j ) = iseed( j )
577  20 CONTINUE
578 *
579 * Initialize RESULT
580 *
581  DO 30 j = 1, 13
582  result( j ) = zero
583  30 CONTINUE
584 *
585 * Generate test matrices A and B
586 *
587 * Description of control parameters:
588 *
589 * KZLASS: =1 means w/o rotation, =2 means w/ rotation,
590 * =3 means random.
591 * KATYPE: the "type" to be passed to ZLATM4 for computing A.
592 * KAZERO: the pattern of zeros on the diagonal for A:
593 * =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
594 * =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
595 * =6: ( 0, 1, 0, xxx, 0 ). (xxx means a string of
596 * non-zero entries.)
597 * KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
598 * =2: large, =3: small.
599 * LASIGN: .TRUE. if the diagonal elements of A are to be
600 * multiplied by a random magnitude 1 number.
601 * KBTYPE, KBZERO, KBMAGN, LBSIGN: the same, but for B.
602 * KTRIAN: =0: don't fill in the upper triangle, =1: do.
603 * KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
604 * RMAGN: used to implement KAMAGN and KBMAGN.
605 *
606  IF( mtypes.GT.maxtyp )
607  $ GO TO 110
608  iinfo = 0
609  IF( kclass( jtype ).LT.3 ) THEN
610 *
611 * Generate A (w/o rotation)
612 *
613  IF( abs( katype( jtype ) ).EQ.3 ) THEN
614  in = 2*( ( n-1 ) / 2 ) + 1
615  IF( in.NE.n )
616  $ CALL zlaset( 'Full', n, n, czero, czero, a, lda )
617  ELSE
618  in = n
619  END IF
620  CALL zlatm4( katype( jtype ), in, kz1( kazero( jtype ) ),
621  $ kz2( kazero( jtype ) ), lasign( jtype ),
622  $ rmagn( kamagn( jtype ) ), ulp,
623  $ rmagn( ktrian( jtype )*kamagn( jtype ) ), 2,
624  $ iseed, a, lda )
625  iadd = kadd( kazero( jtype ) )
626  IF( iadd.GT.0 .AND. iadd.LE.n )
627  $ a( iadd, iadd ) = rmagn( kamagn( jtype ) )
628 *
629 * Generate B (w/o rotation)
630 *
631  IF( abs( kbtype( jtype ) ).EQ.3 ) THEN
632  in = 2*( ( n-1 ) / 2 ) + 1
633  IF( in.NE.n )
634  $ CALL zlaset( 'Full', n, n, czero, czero, b, lda )
635  ELSE
636  in = n
637  END IF
638  CALL zlatm4( kbtype( jtype ), in, kz1( kbzero( jtype ) ),
639  $ kz2( kbzero( jtype ) ), lbsign( jtype ),
640  $ rmagn( kbmagn( jtype ) ), one,
641  $ rmagn( ktrian( jtype )*kbmagn( jtype ) ), 2,
642  $ iseed, b, lda )
643  iadd = kadd( kbzero( jtype ) )
644  IF( iadd.NE.0 .AND. iadd.LE.n )
645  $ b( iadd, iadd ) = rmagn( kbmagn( jtype ) )
646 *
647  IF( kclass( jtype ).EQ.2 .AND. n.GT.0 ) THEN
648 *
649 * Include rotations
650 *
651 * Generate Q, Z as Householder transformations times
652 * a diagonal matrix.
653 *
654  DO 50 jc = 1, n - 1
655  DO 40 jr = jc, n
656  q( jr, jc ) = zlarnd( 3, iseed )
657  z( jr, jc ) = zlarnd( 3, iseed )
658  40 CONTINUE
659  CALL zlarfg( n+1-jc, q( jc, jc ), q( jc+1, jc ), 1,
660  $ work( jc ) )
661  work( 2*n+jc ) = sign( one, dble( q( jc, jc ) ) )
662  q( jc, jc ) = cone
663  CALL zlarfg( n+1-jc, z( jc, jc ), z( jc+1, jc ), 1,
664  $ work( n+jc ) )
665  work( 3*n+jc ) = sign( one, dble( z( jc, jc ) ) )
666  z( jc, jc ) = cone
667  50 CONTINUE
668  ctemp = zlarnd( 3, iseed )
669  q( n, n ) = cone
670  work( n ) = czero
671  work( 3*n ) = ctemp / abs( ctemp )
672  ctemp = zlarnd( 3, iseed )
673  z( n, n ) = cone
674  work( 2*n ) = czero
675  work( 4*n ) = ctemp / abs( ctemp )
676 *
677 * Apply the diagonal matrices
678 *
679  DO 70 jc = 1, n
680  DO 60 jr = 1, n
681  a( jr, jc ) = work( 2*n+jr )*
682  $ dconjg( work( 3*n+jc ) )*
683  $ a( jr, jc )
684  b( jr, jc ) = work( 2*n+jr )*
685  $ dconjg( work( 3*n+jc ) )*
686  $ b( jr, jc )
687  60 CONTINUE
688  70 CONTINUE
689  CALL zunm2r( 'L', 'N', n, n, n-1, q, ldq, work, a,
690  $ lda, work( 2*n+1 ), iinfo )
691  IF( iinfo.NE.0 )
692  $ GO TO 100
693  CALL zunm2r( 'R', 'C', n, n, n-1, z, ldq, work( n+1 ),
694  $ a, lda, work( 2*n+1 ), iinfo )
695  IF( iinfo.NE.0 )
696  $ GO TO 100
697  CALL zunm2r( 'L', 'N', n, n, n-1, q, ldq, work, b,
698  $ lda, work( 2*n+1 ), iinfo )
699  IF( iinfo.NE.0 )
700  $ GO TO 100
701  CALL zunm2r( 'R', 'C', n, n, n-1, z, ldq, work( n+1 ),
702  $ b, lda, work( 2*n+1 ), iinfo )
703  IF( iinfo.NE.0 )
704  $ GO TO 100
705  END IF
706  ELSE
707 *
708 * Random matrices
709 *
710  DO 90 jc = 1, n
711  DO 80 jr = 1, n
712  a( jr, jc ) = rmagn( kamagn( jtype ) )*
713  $ zlarnd( 4, iseed )
714  b( jr, jc ) = rmagn( kbmagn( jtype ) )*
715  $ zlarnd( 4, iseed )
716  80 CONTINUE
717  90 CONTINUE
718  END IF
719 *
720  100 CONTINUE
721 *
722  IF( iinfo.NE.0 ) THEN
723  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
724  $ ioldsd
725  info = abs( iinfo )
726  RETURN
727  END IF
728 *
729  110 CONTINUE
730 *
731  DO 120 i = 1, 13
732  result( i ) = -one
733  120 CONTINUE
734 *
735 * Test with and without sorting of eigenvalues
736 *
737  DO 150 isort = 0, 1
738  IF( isort.EQ.0 ) THEN
739  sort = 'N'
740  rsub = 0
741  ELSE
742  sort = 'S'
743  rsub = 5
744  END IF
745 *
746 * Call ZGGES to compute H, T, Q, Z, alpha, and beta.
747 *
748  CALL zlacpy( 'Full', n, n, a, lda, s, lda )
749  CALL zlacpy( 'Full', n, n, b, lda, t, lda )
750  ntest = 1 + rsub + isort
751  result( 1+rsub+isort ) = ulpinv
752  CALL zgges( 'V', 'V', sort, zlctes, n, s, lda, t, lda,
753  $ sdim, alpha, beta, q, ldq, z, ldq, work,
754  $ lwork, rwork, bwork, iinfo )
755  IF( iinfo.NE.0 .AND. iinfo.NE.n+2 ) THEN
756  result( 1+rsub+isort ) = ulpinv
757  WRITE( nounit, fmt = 9999 )'ZGGES', iinfo, n, jtype,
758  $ ioldsd
759  info = abs( iinfo )
760  GO TO 160
761  END IF
762 *
763  ntest = 4 + rsub
764 *
765 * Do tests 1--4 (or tests 7--9 when reordering )
766 *
767  IF( isort.EQ.0 ) THEN
768  CALL zget51( 1, n, a, lda, s, lda, q, ldq, z, ldq,
769  $ work, rwork, result( 1 ) )
770  CALL zget51( 1, n, b, lda, t, lda, q, ldq, z, ldq,
771  $ work, rwork, result( 2 ) )
772  ELSE
773  CALL zget54( n, a, lda, b, lda, s, lda, t, lda, q,
774  $ ldq, z, ldq, work, result( 2+rsub ) )
775  END IF
776 *
777  CALL zget51( 3, n, b, lda, t, lda, q, ldq, q, ldq, work,
778  $ rwork, result( 3+rsub ) )
779  CALL zget51( 3, n, b, lda, t, lda, z, ldq, z, ldq, work,
780  $ rwork, result( 4+rsub ) )
781 *
782 * Do test 5 and 6 (or Tests 10 and 11 when reordering):
783 * check Schur form of A and compare eigenvalues with
784 * diagonals.
785 *
786  ntest = 6 + rsub
787  temp1 = zero
788 *
789  DO 130 j = 1, n
790  ilabad = .false.
791  temp2 = ( abs1( alpha( j )-s( j, j ) ) /
792  $ max( safmin, abs1( alpha( j ) ), abs1( s( j,
793  $ j ) ) )+abs1( beta( j )-t( j, j ) ) /
794  $ max( safmin, abs1( beta( j ) ), abs1( t( j,
795  $ j ) ) ) ) / ulp
796 *
797  IF( j.LT.n ) THEN
798  IF( s( j+1, j ).NE.zero ) THEN
799  ilabad = .true.
800  result( 5+rsub ) = ulpinv
801  END IF
802  END IF
803  IF( j.GT.1 ) THEN
804  IF( s( j, j-1 ).NE.zero ) THEN
805  ilabad = .true.
806  result( 5+rsub ) = ulpinv
807  END IF
808  END IF
809  temp1 = max( temp1, temp2 )
810  IF( ilabad ) THEN
811  WRITE( nounit, fmt = 9998 )j, n, jtype, ioldsd
812  END IF
813  130 CONTINUE
814  result( 6+rsub ) = temp1
815 *
816  IF( isort.GE.1 ) THEN
817 *
818 * Do test 12
819 *
820  ntest = 12
821  result( 12 ) = zero
822  knteig = 0
823  DO 140 i = 1, n
824  IF( zlctes( alpha( i ), beta( i ) ) )
825  $ knteig = knteig + 1
826  140 CONTINUE
827  IF( sdim.NE.knteig )
828  $ result( 13 ) = ulpinv
829  END IF
830 *
831  150 CONTINUE
832 *
833 * End of Loop -- Check for RESULT(j) > THRESH
834 *
835  160 CONTINUE
836 *
837  ntestt = ntestt + ntest
838 *
839 * Print out tests which fail.
840 *
841  DO 170 jr = 1, ntest
842  IF( result( jr ).GE.thresh ) THEN
843 *
844 * If this is the first test to fail,
845 * print a header to the data file.
846 *
847  IF( nerrs.EQ.0 ) THEN
848  WRITE( nounit, fmt = 9997 )'ZGS'
849 *
850 * Matrix types
851 *
852  WRITE( nounit, fmt = 9996 )
853  WRITE( nounit, fmt = 9995 )
854  WRITE( nounit, fmt = 9994 )'Unitary'
855 *
856 * Tests performed
857 *
858  WRITE( nounit, fmt = 9993 )'unitary', '''',
859  $ 'transpose', ( '''', j = 1, 8 )
860 *
861  END IF
862  nerrs = nerrs + 1
863  IF( result( jr ).LT.10000.0d0 ) THEN
864  WRITE( nounit, fmt = 9992 )n, jtype, ioldsd, jr,
865  $ result( jr )
866  ELSE
867  WRITE( nounit, fmt = 9991 )n, jtype, ioldsd, jr,
868  $ result( jr )
869  END IF
870  END IF
871  170 CONTINUE
872 *
873  180 CONTINUE
874  190 CONTINUE
875 *
876 * Summary
877 *
878  CALL alasvm( 'ZGS', nounit, nerrs, ntestt, 0 )
879 *
880  work( 1 ) = maxwrk
881 *
882  RETURN
883 *
884  9999 FORMAT( ' ZDRGES: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
885  $ i6, ', JTYPE=', i6, ', ISEED=(', 4( i4, ',' ), i5, ')' )
886 *
887  9998 FORMAT( ' ZDRGES: S not in Schur form at eigenvalue ', i6, '.',
888  $ / 9x, 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ),
889  $ i5, ')' )
890 *
891  9997 FORMAT( / 1x, a3, ' -- Complex Generalized Schur from problem ',
892  $ 'driver' )
893 *
894  9996 FORMAT( ' Matrix types (see ZDRGES for details): ' )
895 *
896  9995 FORMAT( ' Special Matrices:', 23x,
897  $ '(J''=transposed Jordan block)',
898  $ / ' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ',
899  $ '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices: ( ',
900  $ 'D=diag(0,1,2,...) )', / ' 7=(D,I) 9=(large*D, small*I',
901  $ ') 11=(large*I, small*D) 13=(large*D, large*I)', /
902  $ ' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ',
903  $ ' 14=(small*D, small*I)', / ' 15=(D, reversed D)' )
904  9994 FORMAT( ' Matrices Rotated by Random ', a, ' Matrices U, V:',
905  $ / ' 16=Transposed Jordan Blocks 19=geometric ',
906  $ 'alpha, beta=0,1', / ' 17=arithm. alpha&beta ',
907  $ ' 20=arithmetic alpha, beta=0,1', / ' 18=clustered ',
908  $ 'alpha, beta=0,1 21=random alpha, beta=0,1',
909  $ / ' Large & Small Matrices:', / ' 22=(large, small) ',
910  $ '23=(small,large) 24=(small,small) 25=(large,large)',
911  $ / ' 26=random O(1) matrices.' )
912 *
913  9993 FORMAT( / ' Tests performed: (S is Schur, T is triangular, ',
914  $ 'Q and Z are ', a, ',', / 19x,
915  $ 'l and r are the appropriate left and right', / 19x,
916  $ 'eigenvectors, resp., a is alpha, b is beta, and', / 19x, a,
917  $ ' means ', a, '.)', / ' Without ordering: ',
918  $ / ' 1 = | A - Q S Z', a,
919  $ ' | / ( |A| n ulp ) 2 = | B - Q T Z', a,
920  $ ' | / ( |B| n ulp )', / ' 3 = | I - QQ', a,
921  $ ' | / ( n ulp ) 4 = | I - ZZ', a,
922  $ ' | / ( n ulp )', / ' 5 = A is in Schur form S',
923  $ / ' 6 = difference between (alpha,beta)',
924  $ ' and diagonals of (S,T)', / ' With ordering: ',
925  $ / ' 7 = | (A,B) - Q (S,T) Z', a, ' | / ( |(A,B)| n ulp )',
926  $ / ' 8 = | I - QQ', a,
927  $ ' | / ( n ulp ) 9 = | I - ZZ', a,
928  $ ' | / ( n ulp )', / ' 10 = A is in Schur form S',
929  $ / ' 11 = difference between (alpha,beta) and diagonals',
930  $ ' of (S,T)', / ' 12 = SDIM is the correct number of ',
931  $ 'selected eigenvalues', / )
932  9992 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
933  $ 4( i4, ',' ), ' result ', i2, ' is', 0p, f8.2 )
934  9991 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
935  $ 4( i4, ',' ), ' result ', i2, ' is', 1p, d10.3 )
936 *
937 * End of ZDRGES
938 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zget54(N, A, LDA, B, LDB, S, LDS, T, LDT, U, LDU, V, LDV, WORK, RESULT)
ZGET54
Definition: zget54.f:158
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:108
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:108
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zgges(JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, WORK, LWORK, RWORK, BWORK, INFO)
ZGGES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition: zgges.f:272
subroutine zunm2r(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
ZUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...
Definition: zunm2r.f:161
subroutine zlatm4(ITYPE, N, NZ1, NZ2, RSIGN, AMAGN, RCOND, TRIANG, IDIST, ISEED, A, LDA)
ZLATM4
Definition: zlatm4.f:173
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
logical function zlctes(Z, D)
ZLCTES
Definition: zlctes.f:60
subroutine zget51(ITYPE, N, A, LDA, B, LDB, U, LDU, V, LDV, WORK, RWORK, RESULT)
ZGET51
Definition: zget51.f:156
complex *16 function zlarnd(IDIST, ISEED)
ZLARND
Definition: zlarnd.f:77

Here is the call graph for this function:

Here is the caller graph for this function: