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

◆ zchkgg()

subroutine zchkgg ( integer  nsizes,
integer, dimension( * )  nn,
integer  ntypes,
logical, dimension( * )  dotype,
integer, dimension( 4 )  iseed,
double precision  thresh,
logical  tstdif,
double precision  thrshn,
integer  nounit,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( lda, * )  b,
complex*16, dimension( lda, * )  h,
complex*16, dimension( lda, * )  t,
complex*16, dimension( lda, * )  s1,
complex*16, dimension( lda, * )  s2,
complex*16, dimension( lda, * )  p1,
complex*16, dimension( lda, * )  p2,
complex*16, dimension( ldu, * )  u,
integer  ldu,
complex*16, dimension( ldu, * )  v,
complex*16, dimension( ldu, * )  q,
complex*16, dimension( ldu, * )  z,
complex*16, dimension( * )  alpha1,
complex*16, dimension( * )  beta1,
complex*16, dimension( * )  alpha3,
complex*16, dimension( * )  beta3,
complex*16, dimension( ldu, * )  evectl,
complex*16, dimension( ldu, * )  evectr,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
logical, dimension( * )  llwork,
double precision, dimension( 15 )  result,
integer  info 
)

ZCHKGG

Purpose:
 ZCHKGG  checks the nonsymmetric generalized eigenvalue problem
 routines.
                                H          H        H
 ZGGHRD factors A and B as U H V  and U T V , where   means conjugate
 transpose, H is hessenberg, T is triangular and U and V are unitary.

                                 H          H
 ZHGEQZ factors H and T as  Q S Z  and Q P Z , where P and S are upper
 triangular and Q and Z are unitary.  It also computes the generalized
 eigenvalues (alpha(1),beta(1)),...,(alpha(n),beta(n)), where
 alpha(j)=S(j,j) and beta(j)=P(j,j) -- thus, w(j) = alpha(j)/beta(j)
 is a root of the generalized eigenvalue problem

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

 and m(j) = beta(j)/alpha(j) is a root of the essentially equivalent
 problem

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

 ZTGEVC computes the matrix L of left eigenvectors and the matrix R
 of right eigenvectors for the matrix pair ( S, P ).  In the
 description below,  l and r are left and right eigenvectors
 corresponding to the generalized eigenvalues (alpha,beta).

 When ZCHKGG is called, a number of matrix "sizes" ("n's") and a
 number of matrix "types" are specified.  For each size ("n")
 and each type of matrix, one matrix will be generated and used
 to test the nonsymmetric eigenroutines.  For each matrix, 13
 tests will be performed.  The first twelve "test ratios" should be
 small -- O(1).  They will be compared with the threshold THRESH:

                  H
 (1)   | A - U H V  | / ( |A| n ulp )

                  H
 (2)   | B - U T V  | / ( |B| n ulp )

               H
 (3)   | I - UU  | / ( n ulp )

               H
 (4)   | I - VV  | / ( n ulp )

                  H
 (5)   | H - Q S Z  | / ( |H| n ulp )

                  H
 (6)   | T - Q P Z  | / ( |T| n ulp )

               H
 (7)   | I - QQ  | / ( n ulp )

               H
 (8)   | I - ZZ  | / ( n ulp )

 (9)   max over all left eigenvalue/-vector pairs (beta/alpha,l) of
                           H
       | (beta A - alpha B) l | / ( ulp max( |beta A|, |alpha B| ) )

 (10)  max over all left eigenvalue/-vector pairs (beta/alpha,l') of
                           H
       | (beta H - alpha T) l' | / ( ulp max( |beta H|, |alpha T| ) )

       where the eigenvectors l' are the result of passing Q to
       DTGEVC and back transforming (JOB='B').

 (11)  max over all right eigenvalue/-vector pairs (beta/alpha,r) of

       | (beta A - alpha B) r | / ( ulp max( |beta A|, |alpha B| ) )

 (12)  max over all right eigenvalue/-vector pairs (beta/alpha,r') of

       | (beta H - alpha T) r' | / ( ulp max( |beta H|, |alpha T| ) )

       where the eigenvectors r' are the result of passing Z to
       DTGEVC and back transforming (JOB='B').

 The last three test ratios will usually be small, but there is no
 mathematical requirement that they be so.  They are therefore
 compared with THRESH only if TSTDIF is .TRUE.

 (13)  | S(Q,Z computed) - S(Q,Z not computed) | / ( |S| ulp )

 (14)  | P(Q,Z computed) - P(Q,Z not computed) | / ( |P| ulp )

 (15)  max( |alpha(Q,Z computed) - alpha(Q,Z not computed)|/|S| ,
            |beta(Q,Z computed) - beta(Q,Z not computed)|/|P| ) / ulp

 In addition, the normalization of L and R are checked, and compared
 with the threshold THRSHN.

 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 P*D1, P is a random unitary diagonal
                       matrix (i.e., with random magnitude 1 entries
                       on the diagonal), and D1=diag( 0, 1,..., N-1 )
                       (i.e., a diagonal matrix with D1(1,1)=0,
                       D1(2,2)=1, ..., D1(N,N)=N-1.)
 (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=P*diag( 0, 0, 1, ..., N-3, 0 ) and
                        D2=Q*diag( 0, N-3, N-4,..., 1, 0, 0 ), and
                        P and Q are random unitary diagonal matrices.
           t   t
 (16) U ( J , J ) V     where U and V are random unitary matrices.

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

 (18) U ( T1, T2 ) V    diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
                        diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
                        s = machine precision.

 (19) U ( T1, T2 ) V    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) U ( T1, T2 ) V    diag(T1)=( 0, 0, 1, 1, a, ..., a   =s, 0 )
                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )

 (21) U ( T1, T2 ) V    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) U ( big*T1, small*T2 ) V   diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
                                 diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (23) U ( small*T1, big*T2 ) V   diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
                                 diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (24) U ( small*T1, small*T2 ) V diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
                                 diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (25) U ( big*T1, big*T2 ) V     diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
                                 diag(T2) = ( 0, 1, ..., 1, 0, 0 )

 (26) U ( T1, T2 ) V     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,
          ZCHKGG 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, ZCHKGG
          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 ZCHKGG 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]TSTDIF
          TSTDIF is LOGICAL
          Specifies whether test ratios 13-15 will be computed and
          compared with THRESH.
          = .FALSE.: Only test ratios 1-12 will be computed and tested.
                     Ratios 13-15 will be set to zero.
          = .TRUE.:  All the test ratios 1-15 will be computed and
                     tested.
[in]THRSHN
          THRSHN is DOUBLE PRECISION
          Threshold for reporting eigenvector normalization error.
          If the normalization of any eigenvector differs from 1 by
          more than THRSHN*ulp, then a special error message will be
          printed.  (This is handled separately from the other tests,
          since only a compiler or programming error should cause an
          error message, at least if THRSHN is at least 5--10.)
[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, H, T, S1, P1, S2, and P2.
          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]H
          H is COMPLEX*16 array, dimension (LDA, max(NN))
          The upper Hessenberg matrix computed from A by ZGGHRD.
[out]T
          T is COMPLEX*16 array, dimension (LDA, max(NN))
          The upper triangular matrix computed from B by ZGGHRD.
[out]S1
          S1 is COMPLEX*16 array, dimension (LDA, max(NN))
          The Schur (upper triangular) matrix computed from H by ZHGEQZ
          when Q and Z are also computed.
[out]S2
          S2 is COMPLEX*16 array, dimension (LDA, max(NN))
          The Schur (upper triangular) matrix computed from H by ZHGEQZ
          when Q and Z are not computed.
[out]P1
          P1 is COMPLEX*16 array, dimension (LDA, max(NN))
          The upper triangular matrix computed from T by ZHGEQZ
          when Q and Z are also computed.
[out]P2
          P2 is COMPLEX*16 array, dimension (LDA, max(NN))
          The upper triangular matrix computed from T by ZHGEQZ
          when Q and Z are not computed.
[out]U
          U is COMPLEX*16 array, dimension (LDU, max(NN))
          The (left) unitary matrix computed by ZGGHRD.
[in]LDU
          LDU is INTEGER
          The leading dimension of U, V, Q, Z, EVECTL, and EVEZTR.  It
          must be at least 1 and at least max( NN ).
[out]V
          V is COMPLEX*16 array, dimension (LDU, max(NN))
          The (right) unitary matrix computed by ZGGHRD.
[out]Q
          Q is COMPLEX*16 array, dimension (LDU, max(NN))
          The (left) unitary matrix computed by ZHGEQZ.
[out]Z
          Z is COMPLEX*16 array, dimension (LDU, max(NN))
          The (left) unitary matrix computed by ZHGEQZ.
[out]ALPHA1
          ALPHA1 is COMPLEX*16 array, dimension (max(NN))
[out]BETA1
          BETA1 is COMPLEX*16 array, dimension (max(NN))
          The generalized eigenvalues of (A,B) computed by ZHGEQZ
          when Q, Z, and the full Schur matrices are computed.
[out]ALPHA3
          ALPHA3 is COMPLEX*16 array, dimension (max(NN))
[out]BETA3
          BETA3 is COMPLEX*16 array, dimension (max(NN))
          The generalized eigenvalues of (A,B) computed by ZHGEQZ
          when neither Q, Z, nor the Schur matrices are computed.
[out]EVECTL
          EVECTL is COMPLEX*16 array, dimension (LDU, max(NN))
          The (lower triangular) left eigenvector matrix for the
          matrices in S1 and P1.
[out]EVECTR
          EVECTR is COMPLEX*16 array, dimension (LDU, max(NN))
          The (upper triangular) right eigenvector matrix for the
          matrices in S1 and P1.
[out]WORK
          WORK is COMPLEX*16 array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          max( 4*N, 2 * N**2, 1 ), for all N=NN(j).
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (2*max(NN))
[out]LLWORK
          LLWORK is LOGICAL array, dimension (max(NN))
[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]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.

Definition at line 498 of file zchkgg.f.

503*
504* -- LAPACK test routine --
505* -- LAPACK is a software package provided by Univ. of Tennessee, --
506* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
507*
508* .. Scalar Arguments ..
509 LOGICAL TSTDIF
510 INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES
511 DOUBLE PRECISION THRESH, THRSHN
512* ..
513* .. Array Arguments ..
514 LOGICAL DOTYPE( * ), LLWORK( * )
515 INTEGER ISEED( 4 ), NN( * )
516 DOUBLE PRECISION RESULT( 15 ), RWORK( * )
517 COMPLEX*16 A( LDA, * ), ALPHA1( * ), ALPHA3( * ),
518 $ B( LDA, * ), BETA1( * ), BETA3( * ),
519 $ EVECTL( LDU, * ), EVECTR( LDU, * ),
520 $ H( LDA, * ), P1( LDA, * ), P2( LDA, * ),
521 $ Q( LDU, * ), S1( LDA, * ), S2( LDA, * ),
522 $ T( LDA, * ), U( LDU, * ), V( LDU, * ),
523 $ WORK( * ), Z( LDU, * )
524* ..
525*
526* =====================================================================
527*
528* .. Parameters ..
529 DOUBLE PRECISION ZERO, ONE
530 parameter( zero = 0.0d+0, one = 1.0d+0 )
531 COMPLEX*16 CZERO, CONE
532 parameter( czero = ( 0.0d+0, 0.0d+0 ),
533 $ cone = ( 1.0d+0, 0.0d+0 ) )
534 INTEGER MAXTYP
535 parameter( maxtyp = 26 )
536* ..
537* .. Local Scalars ..
538 LOGICAL BADNN
539 INTEGER I1, IADD, IINFO, IN, J, JC, JR, JSIZE, JTYPE,
540 $ LWKOPT, MTYPES, N, N1, NERRS, NMATS, NMAX,
541 $ NTEST, NTESTT
542 DOUBLE PRECISION ANORM, BNORM, SAFMAX, SAFMIN, TEMP1, TEMP2,
543 $ ULP, ULPINV
544 COMPLEX*16 CTEMP
545* ..
546* .. Local Arrays ..
547 LOGICAL LASIGN( MAXTYP ), LBSIGN( MAXTYP )
548 INTEGER IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ),
549 $ KATYPE( MAXTYP ), KAZERO( MAXTYP ),
550 $ KBMAGN( MAXTYP ), KBTYPE( MAXTYP ),
551 $ KBZERO( MAXTYP ), KCLASS( MAXTYP ),
552 $ KTRIAN( MAXTYP ), KZ1( 6 ), KZ2( 6 )
553 DOUBLE PRECISION DUMMA( 4 ), RMAGN( 0: 3 )
554 COMPLEX*16 CDUMMA( 4 )
555* ..
556* .. External Functions ..
557 DOUBLE PRECISION DLAMCH, ZLANGE
558 COMPLEX*16 ZLARND
559 EXTERNAL dlamch, zlange, zlarnd
560* ..
561* .. External Subroutines ..
562 EXTERNAL dlasum, xerbla, zgeqr2, zget51, zget52, zgghrd,
564 $ zunm2r
565* ..
566* .. Intrinsic Functions ..
567 INTRINSIC abs, dble, dconjg, max, min, sign
568* ..
569* .. Data statements ..
570 DATA kclass / 15*1, 10*2, 1*3 /
571 DATA kz1 / 0, 1, 2, 1, 3, 3 /
572 DATA kz2 / 0, 0, 1, 2, 1, 1 /
573 DATA kadd / 0, 0, 0, 0, 3, 2 /
574 DATA katype / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
575 $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 /
576 DATA kbtype / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
577 $ 1, 1, -4, 2, -4, 8*8, 0 /
578 DATA kazero / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
579 $ 4*5, 4*3, 1 /
580 DATA kbzero / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
581 $ 4*6, 4*4, 1 /
582 DATA kamagn / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
583 $ 2, 1 /
584 DATA kbmagn / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
585 $ 2, 1 /
586 DATA ktrian / 16*0, 10*1 /
587 DATA lasign / 6*.false., .true., .false., 2*.true.,
588 $ 2*.false., 3*.true., .false., .true.,
589 $ 3*.false., 5*.true., .false. /
590 DATA lbsign / 7*.false., .true., 2*.false.,
591 $ 2*.true., 2*.false., .true., .false., .true.,
592 $ 9*.false. /
593* ..
594* .. Executable Statements ..
595*
596* Check for errors
597*
598 info = 0
599*
600 badnn = .false.
601 nmax = 1
602 DO 10 j = 1, nsizes
603 nmax = max( nmax, nn( j ) )
604 IF( nn( j ).LT.0 )
605 $ badnn = .true.
606 10 CONTINUE
607*
608 lwkopt = max( 2*nmax*nmax, 4*nmax, 1 )
609*
610* Check for errors
611*
612 IF( nsizes.LT.0 ) THEN
613 info = -1
614 ELSE IF( badnn ) THEN
615 info = -2
616 ELSE IF( ntypes.LT.0 ) THEN
617 info = -3
618 ELSE IF( thresh.LT.zero ) THEN
619 info = -6
620 ELSE IF( lda.LE.1 .OR. lda.LT.nmax ) THEN
621 info = -10
622 ELSE IF( ldu.LE.1 .OR. ldu.LT.nmax ) THEN
623 info = -19
624 ELSE IF( lwkopt.GT.lwork ) THEN
625 info = -30
626 END IF
627*
628 IF( info.NE.0 ) THEN
629 CALL xerbla( 'ZCHKGG', -info )
630 RETURN
631 END IF
632*
633* Quick return if possible
634*
635 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
636 $ RETURN
637*
638 safmin = dlamch( 'Safe minimum' )
639 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
640 safmin = safmin / ulp
641 safmax = one / safmin
642 ulpinv = one / ulp
643*
644* The values RMAGN(2:3) depend on N, see below.
645*
646 rmagn( 0 ) = zero
647 rmagn( 1 ) = one
648*
649* Loop over sizes, types
650*
651 ntestt = 0
652 nerrs = 0
653 nmats = 0
654*
655 DO 240 jsize = 1, nsizes
656 n = nn( jsize )
657 n1 = max( 1, n )
658 rmagn( 2 ) = safmax*ulp / dble( n1 )
659 rmagn( 3 ) = safmin*ulpinv*n1
660*
661 IF( nsizes.NE.1 ) THEN
662 mtypes = min( maxtyp, ntypes )
663 ELSE
664 mtypes = min( maxtyp+1, ntypes )
665 END IF
666*
667 DO 230 jtype = 1, mtypes
668 IF( .NOT.dotype( jtype ) )
669 $ GO TO 230
670 nmats = nmats + 1
671 ntest = 0
672*
673* Save ISEED in case of an error.
674*
675 DO 20 j = 1, 4
676 ioldsd( j ) = iseed( j )
677 20 CONTINUE
678*
679* Initialize RESULT
680*
681 DO 30 j = 1, 15
682 result( j ) = zero
683 30 CONTINUE
684*
685* Compute A and B
686*
687* Description of control parameters:
688*
689* KZLASS: =1 means w/o rotation, =2 means w/ rotation,
690* =3 means random.
691* KATYPE: the "type" to be passed to ZLATM4 for computing A.
692* KAZERO: the pattern of zeros on the diagonal for A:
693* =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
694* =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
695* =6: ( 0, 1, 0, xxx, 0 ). (xxx means a string of
696* non-zero entries.)
697* KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
698* =2: large, =3: small.
699* LASIGN: .TRUE. if the diagonal elements of A are to be
700* multiplied by a random magnitude 1 number.
701* KBTYPE, KBZERO, KBMAGN, LBSIGN: the same, but for B.
702* KTRIAN: =0: don't fill in the upper triangle, =1: do.
703* KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
704* RMAGN: used to implement KAMAGN and KBMAGN.
705*
706 IF( mtypes.GT.maxtyp )
707 $ GO TO 110
708 iinfo = 0
709 IF( kclass( jtype ).LT.3 ) THEN
710*
711* Generate A (w/o rotation)
712*
713 IF( abs( katype( jtype ) ).EQ.3 ) THEN
714 in = 2*( ( n-1 ) / 2 ) + 1
715 IF( in.NE.n )
716 $ CALL zlaset( 'Full', n, n, czero, czero, a, lda )
717 ELSE
718 in = n
719 END IF
720 CALL zlatm4( katype( jtype ), in, kz1( kazero( jtype ) ),
721 $ kz2( kazero( jtype ) ), lasign( jtype ),
722 $ rmagn( kamagn( jtype ) ), ulp,
723 $ rmagn( ktrian( jtype )*kamagn( jtype ) ), 4,
724 $ iseed, a, lda )
725 iadd = kadd( kazero( jtype ) )
726 IF( iadd.GT.0 .AND. iadd.LE.n )
727 $ a( iadd, iadd ) = rmagn( kamagn( jtype ) )
728*
729* Generate B (w/o rotation)
730*
731 IF( abs( kbtype( jtype ) ).EQ.3 ) THEN
732 in = 2*( ( n-1 ) / 2 ) + 1
733 IF( in.NE.n )
734 $ CALL zlaset( 'Full', n, n, czero, czero, b, lda )
735 ELSE
736 in = n
737 END IF
738 CALL zlatm4( kbtype( jtype ), in, kz1( kbzero( jtype ) ),
739 $ kz2( kbzero( jtype ) ), lbsign( jtype ),
740 $ rmagn( kbmagn( jtype ) ), one,
741 $ rmagn( ktrian( jtype )*kbmagn( jtype ) ), 4,
742 $ iseed, b, lda )
743 iadd = kadd( kbzero( jtype ) )
744 IF( iadd.NE.0 )
745 $ b( iadd, iadd ) = rmagn( kbmagn( jtype ) )
746*
747 IF( kclass( jtype ).EQ.2 .AND. n.GT.0 ) THEN
748*
749* Include rotations
750*
751* Generate U, V as Householder transformations times a
752* diagonal matrix. (Note that ZLARFG makes U(j,j) and
753* V(j,j) real.)
754*
755 DO 50 jc = 1, n - 1
756 DO 40 jr = jc, n
757 u( jr, jc ) = zlarnd( 3, iseed )
758 v( jr, jc ) = zlarnd( 3, iseed )
759 40 CONTINUE
760 CALL zlarfg( n+1-jc, u( jc, jc ), u( jc+1, jc ), 1,
761 $ work( jc ) )
762 work( 2*n+jc ) = sign( one, dble( u( jc, jc ) ) )
763 u( jc, jc ) = cone
764 CALL zlarfg( n+1-jc, v( jc, jc ), v( jc+1, jc ), 1,
765 $ work( n+jc ) )
766 work( 3*n+jc ) = sign( one, dble( v( jc, jc ) ) )
767 v( jc, jc ) = cone
768 50 CONTINUE
769 ctemp = zlarnd( 3, iseed )
770 u( n, n ) = cone
771 work( n ) = czero
772 work( 3*n ) = ctemp / abs( ctemp )
773 ctemp = zlarnd( 3, iseed )
774 v( n, n ) = cone
775 work( 2*n ) = czero
776 work( 4*n ) = ctemp / abs( ctemp )
777*
778* Apply the diagonal matrices
779*
780 DO 70 jc = 1, n
781 DO 60 jr = 1, n
782 a( jr, jc ) = work( 2*n+jr )*
783 $ dconjg( work( 3*n+jc ) )*
784 $ a( jr, jc )
785 b( jr, jc ) = work( 2*n+jr )*
786 $ dconjg( work( 3*n+jc ) )*
787 $ b( jr, jc )
788 60 CONTINUE
789 70 CONTINUE
790 CALL zunm2r( 'L', 'N', n, n, n-1, u, ldu, work, a,
791 $ lda, work( 2*n+1 ), iinfo )
792 IF( iinfo.NE.0 )
793 $ GO TO 100
794 CALL zunm2r( 'R', 'C', n, n, n-1, v, ldu, work( n+1 ),
795 $ a, lda, work( 2*n+1 ), iinfo )
796 IF( iinfo.NE.0 )
797 $ GO TO 100
798 CALL zunm2r( 'L', 'N', n, n, n-1, u, ldu, work, b,
799 $ lda, work( 2*n+1 ), iinfo )
800 IF( iinfo.NE.0 )
801 $ GO TO 100
802 CALL zunm2r( 'R', 'C', n, n, n-1, v, ldu, work( n+1 ),
803 $ b, lda, work( 2*n+1 ), iinfo )
804 IF( iinfo.NE.0 )
805 $ GO TO 100
806 END IF
807 ELSE
808*
809* Random matrices
810*
811 DO 90 jc = 1, n
812 DO 80 jr = 1, n
813 a( jr, jc ) = rmagn( kamagn( jtype ) )*
814 $ zlarnd( 4, iseed )
815 b( jr, jc ) = rmagn( kbmagn( jtype ) )*
816 $ zlarnd( 4, iseed )
817 80 CONTINUE
818 90 CONTINUE
819 END IF
820*
821 anorm = zlange( '1', n, n, a, lda, rwork )
822 bnorm = zlange( '1', n, n, b, lda, rwork )
823*
824 100 CONTINUE
825*
826 IF( iinfo.NE.0 ) THEN
827 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
828 $ ioldsd
829 info = abs( iinfo )
830 RETURN
831 END IF
832*
833 110 CONTINUE
834*
835* Call ZGEQR2, ZUNM2R, and ZGGHRD to compute H, T, U, and V
836*
837 CALL zlacpy( ' ', n, n, a, lda, h, lda )
838 CALL zlacpy( ' ', n, n, b, lda, t, lda )
839 ntest = 1
840 result( 1 ) = ulpinv
841*
842 CALL zgeqr2( n, n, t, lda, work, work( n+1 ), iinfo )
843 IF( iinfo.NE.0 ) THEN
844 WRITE( nounit, fmt = 9999 )'ZGEQR2', iinfo, n, jtype,
845 $ ioldsd
846 info = abs( iinfo )
847 GO TO 210
848 END IF
849*
850 CALL zunm2r( 'L', 'C', n, n, n, t, lda, work, h, lda,
851 $ work( n+1 ), iinfo )
852 IF( iinfo.NE.0 ) THEN
853 WRITE( nounit, fmt = 9999 )'ZUNM2R', iinfo, n, jtype,
854 $ ioldsd
855 info = abs( iinfo )
856 GO TO 210
857 END IF
858*
859 CALL zlaset( 'Full', n, n, czero, cone, u, ldu )
860 CALL zunm2r( 'R', 'N', n, n, n, t, lda, work, u, ldu,
861 $ work( n+1 ), iinfo )
862 IF( iinfo.NE.0 ) THEN
863 WRITE( nounit, fmt = 9999 )'ZUNM2R', iinfo, n, jtype,
864 $ ioldsd
865 info = abs( iinfo )
866 GO TO 210
867 END IF
868*
869 CALL zgghrd( 'V', 'I', n, 1, n, h, lda, t, lda, u, ldu, v,
870 $ ldu, iinfo )
871 IF( iinfo.NE.0 ) THEN
872 WRITE( nounit, fmt = 9999 )'ZGGHRD', iinfo, n, jtype,
873 $ ioldsd
874 info = abs( iinfo )
875 GO TO 210
876 END IF
877 ntest = 4
878*
879* Do tests 1--4
880*
881 CALL zget51( 1, n, a, lda, h, lda, u, ldu, v, ldu, work,
882 $ rwork, result( 1 ) )
883 CALL zget51( 1, n, b, lda, t, lda, u, ldu, v, ldu, work,
884 $ rwork, result( 2 ) )
885 CALL zget51( 3, n, b, lda, t, lda, u, ldu, u, ldu, work,
886 $ rwork, result( 3 ) )
887 CALL zget51( 3, n, b, lda, t, lda, v, ldu, v, ldu, work,
888 $ rwork, result( 4 ) )
889*
890* Call ZHGEQZ to compute S1, P1, S2, P2, Q, and Z, do tests.
891*
892* Compute T1 and UZ
893*
894* Eigenvalues only
895*
896 CALL zlacpy( ' ', n, n, h, lda, s2, lda )
897 CALL zlacpy( ' ', n, n, t, lda, p2, lda )
898 ntest = 5
899 result( 5 ) = ulpinv
900*
901 CALL zhgeqz( 'E', 'N', 'N', n, 1, n, s2, lda, p2, lda,
902 $ alpha3, beta3, q, ldu, z, ldu, work, lwork,
903 $ rwork, iinfo )
904 IF( iinfo.NE.0 ) THEN
905 WRITE( nounit, fmt = 9999 )'ZHGEQZ(E)', iinfo, n, jtype,
906 $ ioldsd
907 info = abs( iinfo )
908 GO TO 210
909 END IF
910*
911* Eigenvalues and Full Schur Form
912*
913 CALL zlacpy( ' ', n, n, h, lda, s2, lda )
914 CALL zlacpy( ' ', n, n, t, lda, p2, lda )
915*
916 CALL zhgeqz( 'S', 'N', 'N', n, 1, n, s2, lda, p2, lda,
917 $ alpha1, beta1, q, ldu, z, ldu, work, lwork,
918 $ rwork, iinfo )
919 IF( iinfo.NE.0 ) THEN
920 WRITE( nounit, fmt = 9999 )'ZHGEQZ(S)', iinfo, n, jtype,
921 $ ioldsd
922 info = abs( iinfo )
923 GO TO 210
924 END IF
925*
926* Eigenvalues, Schur Form, and Schur Vectors
927*
928 CALL zlacpy( ' ', n, n, h, lda, s1, lda )
929 CALL zlacpy( ' ', n, n, t, lda, p1, lda )
930*
931 CALL zhgeqz( 'S', 'I', 'I', n, 1, n, s1, lda, p1, lda,
932 $ alpha1, beta1, q, ldu, z, ldu, work, lwork,
933 $ rwork, iinfo )
934 IF( iinfo.NE.0 ) THEN
935 WRITE( nounit, fmt = 9999 )'ZHGEQZ(V)', iinfo, n, jtype,
936 $ ioldsd
937 info = abs( iinfo )
938 GO TO 210
939 END IF
940*
941 ntest = 8
942*
943* Do Tests 5--8
944*
945 CALL zget51( 1, n, h, lda, s1, lda, q, ldu, z, ldu, work,
946 $ rwork, result( 5 ) )
947 CALL zget51( 1, n, t, lda, p1, lda, q, ldu, z, ldu, work,
948 $ rwork, result( 6 ) )
949 CALL zget51( 3, n, t, lda, p1, lda, q, ldu, q, ldu, work,
950 $ rwork, result( 7 ) )
951 CALL zget51( 3, n, t, lda, p1, lda, z, ldu, z, ldu, work,
952 $ rwork, result( 8 ) )
953*
954* Compute the Left and Right Eigenvectors of (S1,P1)
955*
956* 9: Compute the left eigenvector Matrix without
957* back transforming:
958*
959 ntest = 9
960 result( 9 ) = ulpinv
961*
962* To test "SELECT" option, compute half of the eigenvectors
963* in one call, and half in another
964*
965 i1 = n / 2
966 DO 120 j = 1, i1
967 llwork( j ) = .true.
968 120 CONTINUE
969 DO 130 j = i1 + 1, n
970 llwork( j ) = .false.
971 130 CONTINUE
972*
973 CALL ztgevc( 'L', 'S', llwork, n, s1, lda, p1, lda, evectl,
974 $ ldu, cdumma, ldu, n, in, work, rwork, iinfo )
975 IF( iinfo.NE.0 ) THEN
976 WRITE( nounit, fmt = 9999 )'ZTGEVC(L,S1)', iinfo, n,
977 $ jtype, ioldsd
978 info = abs( iinfo )
979 GO TO 210
980 END IF
981*
982 i1 = in
983 DO 140 j = 1, i1
984 llwork( j ) = .false.
985 140 CONTINUE
986 DO 150 j = i1 + 1, n
987 llwork( j ) = .true.
988 150 CONTINUE
989*
990 CALL ztgevc( 'L', 'S', llwork, n, s1, lda, p1, lda,
991 $ evectl( 1, i1+1 ), ldu, cdumma, ldu, n, in,
992 $ work, rwork, iinfo )
993 IF( iinfo.NE.0 ) THEN
994 WRITE( nounit, fmt = 9999 )'ZTGEVC(L,S2)', iinfo, n,
995 $ jtype, ioldsd
996 info = abs( iinfo )
997 GO TO 210
998 END IF
999*
1000 CALL zget52( .true., n, s1, lda, p1, lda, evectl, ldu,
1001 $ alpha1, beta1, work, rwork, dumma( 1 ) )
1002 result( 9 ) = dumma( 1 )
1003 IF( dumma( 2 ).GT.thrshn ) THEN
1004 WRITE( nounit, fmt = 9998 )'Left', 'ZTGEVC(HOWMNY=S)',
1005 $ dumma( 2 ), n, jtype, ioldsd
1006 END IF
1007*
1008* 10: Compute the left eigenvector Matrix with
1009* back transforming:
1010*
1011 ntest = 10
1012 result( 10 ) = ulpinv
1013 CALL zlacpy( 'F', n, n, q, ldu, evectl, ldu )
1014 CALL ztgevc( 'L', 'B', llwork, n, s1, lda, p1, lda, evectl,
1015 $ ldu, cdumma, ldu, n, in, work, rwork, iinfo )
1016 IF( iinfo.NE.0 ) THEN
1017 WRITE( nounit, fmt = 9999 )'ZTGEVC(L,B)', iinfo, n,
1018 $ jtype, ioldsd
1019 info = abs( iinfo )
1020 GO TO 210
1021 END IF
1022*
1023 CALL zget52( .true., n, h, lda, t, lda, evectl, ldu, alpha1,
1024 $ beta1, work, rwork, dumma( 1 ) )
1025 result( 10 ) = dumma( 1 )
1026 IF( dumma( 2 ).GT.thrshn ) THEN
1027 WRITE( nounit, fmt = 9998 )'Left', 'ZTGEVC(HOWMNY=B)',
1028 $ dumma( 2 ), n, jtype, ioldsd
1029 END IF
1030*
1031* 11: Compute the right eigenvector Matrix without
1032* back transforming:
1033*
1034 ntest = 11
1035 result( 11 ) = ulpinv
1036*
1037* To test "SELECT" option, compute half of the eigenvectors
1038* in one call, and half in another
1039*
1040 i1 = n / 2
1041 DO 160 j = 1, i1
1042 llwork( j ) = .true.
1043 160 CONTINUE
1044 DO 170 j = i1 + 1, n
1045 llwork( j ) = .false.
1046 170 CONTINUE
1047*
1048 CALL ztgevc( 'R', 'S', llwork, n, s1, lda, p1, lda, cdumma,
1049 $ ldu, evectr, ldu, n, in, work, rwork, iinfo )
1050 IF( iinfo.NE.0 ) THEN
1051 WRITE( nounit, fmt = 9999 )'ZTGEVC(R,S1)', iinfo, n,
1052 $ jtype, ioldsd
1053 info = abs( iinfo )
1054 GO TO 210
1055 END IF
1056*
1057 i1 = in
1058 DO 180 j = 1, i1
1059 llwork( j ) = .false.
1060 180 CONTINUE
1061 DO 190 j = i1 + 1, n
1062 llwork( j ) = .true.
1063 190 CONTINUE
1064*
1065 CALL ztgevc( 'R', 'S', llwork, n, s1, lda, p1, lda, cdumma,
1066 $ ldu, evectr( 1, i1+1 ), ldu, n, in, work,
1067 $ rwork, iinfo )
1068 IF( iinfo.NE.0 ) THEN
1069 WRITE( nounit, fmt = 9999 )'ZTGEVC(R,S2)', iinfo, n,
1070 $ jtype, ioldsd
1071 info = abs( iinfo )
1072 GO TO 210
1073 END IF
1074*
1075 CALL zget52( .false., n, s1, lda, p1, lda, evectr, ldu,
1076 $ alpha1, beta1, work, rwork, dumma( 1 ) )
1077 result( 11 ) = dumma( 1 )
1078 IF( dumma( 2 ).GT.thresh ) THEN
1079 WRITE( nounit, fmt = 9998 )'Right', 'ZTGEVC(HOWMNY=S)',
1080 $ dumma( 2 ), n, jtype, ioldsd
1081 END IF
1082*
1083* 12: Compute the right eigenvector Matrix with
1084* back transforming:
1085*
1086 ntest = 12
1087 result( 12 ) = ulpinv
1088 CALL zlacpy( 'F', n, n, z, ldu, evectr, ldu )
1089 CALL ztgevc( 'R', 'B', llwork, n, s1, lda, p1, lda, cdumma,
1090 $ ldu, evectr, ldu, n, in, work, rwork, iinfo )
1091 IF( iinfo.NE.0 ) THEN
1092 WRITE( nounit, fmt = 9999 )'ZTGEVC(R,B)', iinfo, n,
1093 $ jtype, ioldsd
1094 info = abs( iinfo )
1095 GO TO 210
1096 END IF
1097*
1098 CALL zget52( .false., n, h, lda, t, lda, evectr, ldu,
1099 $ alpha1, beta1, work, rwork, dumma( 1 ) )
1100 result( 12 ) = dumma( 1 )
1101 IF( dumma( 2 ).GT.thresh ) THEN
1102 WRITE( nounit, fmt = 9998 )'Right', 'ZTGEVC(HOWMNY=B)',
1103 $ dumma( 2 ), n, jtype, ioldsd
1104 END IF
1105*
1106* Tests 13--15 are done only on request
1107*
1108 IF( tstdif ) THEN
1109*
1110* Do Tests 13--14
1111*
1112 CALL zget51( 2, n, s1, lda, s2, lda, q, ldu, z, ldu,
1113 $ work, rwork, result( 13 ) )
1114 CALL zget51( 2, n, p1, lda, p2, lda, q, ldu, z, ldu,
1115 $ work, rwork, result( 14 ) )
1116*
1117* Do Test 15
1118*
1119 temp1 = zero
1120 temp2 = zero
1121 DO 200 j = 1, n
1122 temp1 = max( temp1, abs( alpha1( j )-alpha3( j ) ) )
1123 temp2 = max( temp2, abs( beta1( j )-beta3( j ) ) )
1124 200 CONTINUE
1125*
1126 temp1 = temp1 / max( safmin, ulp*max( temp1, anorm ) )
1127 temp2 = temp2 / max( safmin, ulp*max( temp2, bnorm ) )
1128 result( 15 ) = max( temp1, temp2 )
1129 ntest = 15
1130 ELSE
1131 result( 13 ) = zero
1132 result( 14 ) = zero
1133 result( 15 ) = zero
1134 ntest = 12
1135 END IF
1136*
1137* End of Loop -- Check for RESULT(j) > THRESH
1138*
1139 210 CONTINUE
1140*
1141 ntestt = ntestt + ntest
1142*
1143* Print out tests which fail.
1144*
1145 DO 220 jr = 1, ntest
1146 IF( result( jr ).GE.thresh ) THEN
1147*
1148* If this is the first test to fail,
1149* print a header to the data file.
1150*
1151 IF( nerrs.EQ.0 ) THEN
1152 WRITE( nounit, fmt = 9997 )'ZGG'
1153*
1154* Matrix types
1155*
1156 WRITE( nounit, fmt = 9996 )
1157 WRITE( nounit, fmt = 9995 )
1158 WRITE( nounit, fmt = 9994 )'Unitary'
1159*
1160* Tests performed
1161*
1162 WRITE( nounit, fmt = 9993 )'unitary', '*',
1163 $ 'conjugate transpose', ( '*', j = 1, 10 )
1164*
1165 END IF
1166 nerrs = nerrs + 1
1167 IF( result( jr ).LT.10000.0d0 ) THEN
1168 WRITE( nounit, fmt = 9992 )n, jtype, ioldsd, jr,
1169 $ result( jr )
1170 ELSE
1171 WRITE( nounit, fmt = 9991 )n, jtype, ioldsd, jr,
1172 $ result( jr )
1173 END IF
1174 END IF
1175 220 CONTINUE
1176*
1177 230 CONTINUE
1178 240 CONTINUE
1179*
1180* Summary
1181*
1182 CALL dlasum( 'ZGG', nounit, nerrs, ntestt )
1183 RETURN
1184*
1185 9999 FORMAT( ' ZCHKGG: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1186 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1187*
1188 9998 FORMAT( ' ZCHKGG: ', a, ' Eigenvectors from ', a, ' incorrectly ',
1189 $ 'normalized.', / ' Bits of error=', 0p, g10.3, ',', 9x,
1190 $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
1191 $ ')' )
1192*
1193 9997 FORMAT( 1x, a3, ' -- Complex Generalized eigenvalue problem' )
1194*
1195 9996 FORMAT( ' Matrix types (see ZCHKGG for details): ' )
1196*
1197 9995 FORMAT( ' Special Matrices:', 23x,
1198 $ '(J''=transposed Jordan block)',
1199 $ / ' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ',
1200 $ '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices: ( ',
1201 $ 'D=diag(0,1,2,...) )', / ' 7=(D,I) 9=(large*D, small*I',
1202 $ ') 11=(large*I, small*D) 13=(large*D, large*I)', /
1203 $ ' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ',
1204 $ ' 14=(small*D, small*I)', / ' 15=(D, reversed D)' )
1205 9994 FORMAT( ' Matrices Rotated by Random ', a, ' Matrices U, V:',
1206 $ / ' 16=Transposed Jordan Blocks 19=geometric ',
1207 $ 'alpha, beta=0,1', / ' 17=arithm. alpha&beta ',
1208 $ ' 20=arithmetic alpha, beta=0,1', / ' 18=clustered ',
1209 $ 'alpha, beta=0,1 21=random alpha, beta=0,1',
1210 $ / ' Large & Small Matrices:', / ' 22=(large, small) ',
1211 $ '23=(small,large) 24=(small,small) 25=(large,large)',
1212 $ / ' 26=random O(1) matrices.' )
1213*
1214 9993 FORMAT( / ' Tests performed: (H is Hessenberg, S is Schur, B, ',
1215 $ 'T, P are triangular,', / 20x, 'U, V, Q, and Z are ', a,
1216 $ ', l and r are the', / 20x,
1217 $ 'appropriate left and right eigenvectors, resp., a is',
1218 $ / 20x, 'alpha, b is beta, and ', a, ' means ', a, '.)',
1219 $ / ' 1 = | A - U H V', a,
1220 $ ' | / ( |A| n ulp ) 2 = | B - U T V', a,
1221 $ ' | / ( |B| n ulp )', / ' 3 = | I - UU', a,
1222 $ ' | / ( n ulp ) 4 = | I - VV', a,
1223 $ ' | / ( n ulp )', / ' 5 = | H - Q S Z', a,
1224 $ ' | / ( |H| n ulp )', 6x, '6 = | T - Q P Z', a,
1225 $ ' | / ( |T| n ulp )', / ' 7 = | I - QQ', a,
1226 $ ' | / ( n ulp ) 8 = | I - ZZ', a,
1227 $ ' | / ( n ulp )', / ' 9 = max | ( b S - a P )', a,
1228 $ ' l | / const. 10 = max | ( b H - a T )', a,
1229 $ ' l | / const.', /
1230 $ ' 11= max | ( b S - a P ) r | / const. 12 = max | ( b H',
1231 $ ' - a T ) r | / const.', / 1x )
1232*
1233 9992 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
1234 $ 4( i4, ',' ), ' result ', i2, ' is', 0p, f8.2 )
1235 9991 FORMAT( ' Matrix order=', i5, ', type=', i2, ', seed=',
1236 $ 4( i4, ',' ), ' result ', i2, ' is', 1p, d10.3 )
1237*
1238* End of ZCHKGG
1239*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlasum(type, iounit, ie, nrun)
DLASUM
Definition dlasum.f:43
subroutine zgeqr2(m, n, a, lda, tau, work, info)
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgeqr2.f:130
subroutine zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:204
subroutine zhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
ZHGEQZ
Definition zhgeqz.f:284
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:106
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 ztgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
ZTGEVC
Definition ztgevc.f:219
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:159
subroutine zget51(itype, n, a, lda, b, ldb, u, ldu, v, ldv, work, rwork, result)
ZGET51
Definition zget51.f:155
subroutine zget52(left, n, a, lda, b, ldb, e, lde, alpha, beta, work, rwork, result)
ZGET52
Definition zget52.f:162
complex *16 function zlarnd(idist, iseed)
ZLARND
Definition zlarnd.f:75
subroutine zlatm4(itype, n, nz1, nz2, rsign, amagn, rcond, triang, idist, iseed, a, lda)
ZLATM4
Definition zlatm4.f:171
Here is the call graph for this function:
Here is the caller graph for this function: