3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
11 * SUBROUTINE DCHKGG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12 * TSTDIF, THRSHN, NOUNIT, A, LDA, B, H, T, S1,
13 * S2, P1, P2, U, LDU, V, Q, Z, ALPHR1, ALPHI1,
14 * BETA1, ALPHR3, ALPHI3, BETA3, EVECTL, EVECTR,
15 * WORK, LWORK, LLWORK, RESULT, INFO )
17 * .. Scalar Arguments ..
19 * INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES
20 * DOUBLE PRECISION THRESH, THRSHN
22 * .. Array Arguments ..
23 * LOGICAL DOTYPE( * ), LLWORK( * )
24 * INTEGER ISEED( 4 ), NN( * )
25 * DOUBLE PRECISION A( LDA, * ), ALPHI1( * ), ALPHI3( * ),
26 * $ ALPHR1( * ), ALPHR3( * ), B( LDA, * ),
27 * $ BETA1( * ), BETA3( * ), EVECTL( LDU, * ),
28 * $ EVECTR( LDU, * ), H( LDA, * ), P1( LDA, * ),
29 * $ P2( LDA, * ), Q( LDU, * ), RESULT( 15 ),
30 * $ S1( LDA, * ), S2( LDA, * ), T( LDA, * ),
31 * $ U( LDU, * ), V( LDU, * ), WORK( * ),
41 *> DCHKGG checks the nonsymmetric generalized eigenvalue problem
44 *> DGGHRD factors A and B as U H V and U T V , where means
45 *> transpose, H is hessenberg, T is triangular and U and V are
48 *> DHGEQZ factors H and T as Q S Z and Q P Z , where P is upper
49 *> triangular, S is in generalized Schur form (block upper triangular,
50 *> with 1x1 and 2x2 blocks on the diagonal, the 2x2 blocks
51 *> corresponding to complex conjugate pairs of generalized
52 *> eigenvalues), and Q and Z are orthogonal. It also computes the
53 *> generalized eigenvalues (alpha(1),beta(1)),...,(alpha(n),beta(n)),
54 *> where alpha(j)=S(j,j) and beta(j)=P(j,j) -- thus,
55 *> w(j) = alpha(j)/beta(j) is a root of the generalized eigenvalue
58 *> det( A - w(j) B ) = 0
60 *> and m(j) = beta(j)/alpha(j) is a root of the essentially equivalent
63 *> det( m(j) A - B ) = 0
65 *> DTGEVC computes the matrix L of left eigenvectors and the matrix R
66 *> of right eigenvectors for the matrix pair ( S, P ). In the
67 *> description below, l and r are left and right eigenvectors
68 *> corresponding to the generalized eigenvalues (alpha,beta).
70 *> When DCHKGG is called, a number of matrix "sizes" ("n's") and a
71 *> number of matrix "types" are specified. For each size ("n")
72 *> and each type of matrix, one matrix will be generated and used
73 *> to test the nonsymmetric eigenroutines. For each matrix, 15
74 *> tests will be performed. The first twelve "test ratios" should be
75 *> small -- O(1). They will be compared with the threshold THRESH:
78 *> (1) | A - U H V | / ( |A| n ulp )
81 *> (2) | B - U T V | / ( |B| n ulp )
84 *> (3) | I - UU | / ( n ulp )
87 *> (4) | I - VV | / ( n ulp )
90 *> (5) | H - Q S Z | / ( |H| n ulp )
93 *> (6) | T - Q P Z | / ( |T| n ulp )
96 *> (7) | I - QQ | / ( n ulp )
99 *> (8) | I - ZZ | / ( n ulp )
101 *> (9) max over all left eigenvalue/-vector pairs (beta/alpha,l) of
103 *> | l**H * (beta S - alpha P) | / ( ulp max( |beta S|, |alpha P| ) )
105 *> (10) max over all left eigenvalue/-vector pairs (beta/alpha,l') of
107 *> | l'**H * (beta H - alpha T) | / ( ulp max( |beta H|, |alpha T| ) )
109 *> where the eigenvectors l' are the result of passing Q to
110 *> DTGEVC and back transforming (HOWMNY='B').
112 *> (11) max over all right eigenvalue/-vector pairs (beta/alpha,r) of
114 *> | (beta S - alpha T) r | / ( ulp max( |beta S|, |alpha T| ) )
116 *> (12) max over all right eigenvalue/-vector pairs (beta/alpha,r') of
118 *> | (beta H - alpha T) r' | / ( ulp max( |beta H|, |alpha T| ) )
120 *> where the eigenvectors r' are the result of passing Z to
121 *> DTGEVC and back transforming (HOWMNY='B').
123 *> The last three test ratios will usually be small, but there is no
124 *> mathematical requirement that they be so. They are therefore
125 *> compared with THRESH only if TSTDIF is .TRUE.
127 *> (13) | S(Q,Z computed) - S(Q,Z not computed) | / ( |S| ulp )
129 *> (14) | P(Q,Z computed) - P(Q,Z not computed) | / ( |P| ulp )
131 *> (15) max( |alpha(Q,Z computed) - alpha(Q,Z not computed)|/|S| ,
132 *> |beta(Q,Z computed) - beta(Q,Z not computed)|/|P| ) / ulp
134 *> In addition, the normalization of L and R are checked, and compared
135 *> with the threshold THRSHN.
140 *> The sizes of the test matrices are specified by an array
141 *> NN(1:NSIZES); the value of each element NN(j) specifies one size.
142 *> The "types" are specified by a logical array DOTYPE( 1:NTYPES ); if
143 *> DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
144 *> Currently, the list of possible types is:
146 *> (1) ( 0, 0 ) (a pair of zero matrices)
148 *> (2) ( I, 0 ) (an identity and a zero matrix)
150 *> (3) ( 0, I ) (an identity and a zero matrix)
152 *> (4) ( I, I ) (a pair of identity matrices)
155 *> (5) ( J , J ) (a pair of transposed Jordan blocks)
158 *> (6) ( X, Y ) where X = ( J 0 ) and Y = ( t )
160 *> and I is a k x k identity and J a (k+1)x(k+1)
161 *> Jordan block; k=(N-1)/2
163 *> (7) ( D, I ) where D is diag( 0, 1,..., N-1 ) (a diagonal
164 *> matrix with those diagonal entries.)
167 *> (9) ( big*D, small*I ) where "big" is near overflow and small=1/big
169 *> (10) ( small*D, big*I )
171 *> (11) ( big*I, small*D )
173 *> (12) ( small*I, big*D )
175 *> (13) ( big*D, big*I )
177 *> (14) ( small*D, small*I )
179 *> (15) ( D1, D2 ) where D1 is diag( 0, 0, 1, ..., N-3, 0 ) and
180 *> D2 is diag( 0, N-3, N-4,..., 1, 0, 0 )
182 *> (16) U ( J , J ) V where U and V are random orthogonal matrices.
184 *> (17) U ( T1, T2 ) V where T1 and T2 are upper triangular matrices
185 *> with random O(1) entries above the diagonal
186 *> and diagonal entries diag(T1) =
187 *> ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
188 *> ( 0, N-3, N-4,..., 1, 0, 0 )
190 *> (18) U ( T1, T2 ) V diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
191 *> diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
192 *> s = machine precision.
194 *> (19) U ( T1, T2 ) V diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 )
195 *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 )
198 *> (20) U ( T1, T2 ) V diag(T1)=( 0, 0, 1, 1, a, ..., a =s, 0 )
199 *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
201 *> (21) U ( T1, T2 ) V diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 )
202 *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
203 *> where r1,..., r(N-4) are random.
205 *> (22) U ( big*T1, small*T2 ) V diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
206 *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
208 *> (23) U ( small*T1, big*T2 ) V diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
209 *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
211 *> (24) U ( small*T1, small*T2 ) V diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
212 *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
214 *> (25) U ( big*T1, big*T2 ) V diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
215 *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
217 *> (26) U ( T1, T2 ) V where T1 and T2 are random upper-triangular
227 *> The number of sizes of matrices to use. If it is zero,
228 *> DCHKGG does nothing. It must be at least zero.
233 *> NN is INTEGER array, dimension (NSIZES)
234 *> An array containing the sizes to be used for the matrices.
235 *> Zero values will be skipped. The values must be at least
242 *> The number of elements in DOTYPE. If it is zero, DCHKGG
243 *> does nothing. It must be at least zero. If it is MAXTYP+1
244 *> and NSIZES is 1, then an additional type, MAXTYP+1 is
245 *> defined, which is to use whatever matrix is in A. This
246 *> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
247 *> DOTYPE(MAXTYP+1) is .TRUE. .
252 *> DOTYPE is LOGICAL array, dimension (NTYPES)
253 *> If DOTYPE(j) is .TRUE., then for each size in NN a
254 *> matrix of that size and of type j will be generated.
255 *> If NTYPES is smaller than the maximum number of types
256 *> defined (PARAMETER MAXTYP), then types NTYPES+1 through
257 *> MAXTYP will not be generated. If NTYPES is larger
258 *> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
262 *> \param[in,out] ISEED
264 *> ISEED is INTEGER array, dimension (4)
265 *> On entry ISEED specifies the seed of the random number
266 *> generator. The array elements should be between 0 and 4095;
267 *> if not they will be reduced mod 4096. Also, ISEED(4) must
268 *> be odd. The random number generator uses a linear
269 *> congruential sequence limited to small integers, and so
270 *> should produce machine independent random numbers. The
271 *> values of ISEED are changed on exit, and can be used in the
272 *> next call to DCHKGG to continue the same random number
278 *> THRESH is DOUBLE PRECISION
279 *> A test will count as "failed" if the "error", computed as
280 *> described above, exceeds THRESH. Note that the error is
281 *> scaled to be O(1), so THRESH should be a reasonably small
282 *> multiple of 1, e.g., 10 or 100. In particular, it should
283 *> not depend on the precision (single vs. double) or the size
284 *> of the matrix. It must be at least zero.
290 *> Specifies whether test ratios 13-15 will be computed and
291 *> compared with THRESH.
292 *> = .FALSE.: Only test ratios 1-12 will be computed and tested.
293 *> Ratios 13-15 will be set to zero.
294 *> = .TRUE.: All the test ratios 1-15 will be computed and
300 *> THRSHN is DOUBLE PRECISION
301 *> Threshold for reporting eigenvector normalization error.
302 *> If the normalization of any eigenvector differs from 1 by
303 *> more than THRSHN*ulp, then a special error message will be
304 *> printed. (This is handled separately from the other tests,
305 *> since only a compiler or programming error should cause an
306 *> error message, at least if THRSHN is at least 5--10.)
312 *> The FORTRAN unit number for printing out error messages
313 *> (e.g., if a routine returns IINFO not equal to 0.)
318 *> A is DOUBLE PRECISION array, dimension
320 *> Used to hold the original A matrix. Used as input only
321 *> if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
322 *> DOTYPE(MAXTYP+1)=.TRUE.
328 *> The leading dimension of A, B, H, T, S1, P1, S2, and P2.
329 *> It must be at least 1 and at least max( NN ).
334 *> B is DOUBLE PRECISION array, dimension
336 *> Used to hold the original B matrix. Used as input only
337 *> if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
338 *> DOTYPE(MAXTYP+1)=.TRUE.
343 *> H is DOUBLE PRECISION array, dimension (LDA, max(NN))
344 *> The upper Hessenberg matrix computed from A by DGGHRD.
349 *> T is DOUBLE PRECISION array, dimension (LDA, max(NN))
350 *> The upper triangular matrix computed from B by DGGHRD.
355 *> S1 is DOUBLE PRECISION array, dimension (LDA, max(NN))
356 *> The Schur (block upper triangular) matrix computed from H by
357 *> DHGEQZ when Q and Z are also computed.
362 *> S2 is DOUBLE PRECISION array, dimension (LDA, max(NN))
363 *> The Schur (block upper triangular) matrix computed from H by
364 *> DHGEQZ when Q and Z are not computed.
369 *> P1 is DOUBLE PRECISION array, dimension (LDA, max(NN))
370 *> The upper triangular matrix computed from T by DHGEQZ
371 *> when Q and Z are also computed.
376 *> P2 is DOUBLE PRECISION array, dimension (LDA, max(NN))
377 *> The upper triangular matrix computed from T by DHGEQZ
378 *> when Q and Z are not computed.
383 *> U is DOUBLE PRECISION array, dimension (LDU, max(NN))
384 *> The (left) orthogonal matrix computed by DGGHRD.
390 *> The leading dimension of U, V, Q, Z, EVECTL, and EVEZTR. It
391 *> must be at least 1 and at least max( NN ).
396 *> V is DOUBLE PRECISION array, dimension (LDU, max(NN))
397 *> The (right) orthogonal matrix computed by DGGHRD.
402 *> Q is DOUBLE PRECISION array, dimension (LDU, max(NN))
403 *> The (left) orthogonal matrix computed by DHGEQZ.
408 *> Z is DOUBLE PRECISION array, dimension (LDU, max(NN))
409 *> The (left) orthogonal matrix computed by DHGEQZ.
412 *> \param[out] ALPHR1
414 *> ALPHR1 is DOUBLE PRECISION array, dimension (max(NN))
417 *> \param[out] ALPHI1
419 *> ALPHI1 is DOUBLE PRECISION array, dimension (max(NN))
424 *> BETA1 is DOUBLE PRECISION array, dimension (max(NN))
426 *> The generalized eigenvalues of (A,B) computed by DHGEQZ
427 *> when Q, Z, and the full Schur matrices are computed.
428 *> On exit, ( ALPHR1(k)+ALPHI1(k)*i ) / BETA1(k) is the k-th
429 *> generalized eigenvalue of the matrices in A and B.
432 *> \param[out] ALPHR3
434 *> ALPHR3 is DOUBLE PRECISION array, dimension (max(NN))
437 *> \param[out] ALPHI3
439 *> ALPHI3 is DOUBLE PRECISION array, dimension (max(NN))
444 *> BETA3 is DOUBLE PRECISION array, dimension (max(NN))
447 *> \param[out] EVECTL
449 *> EVECTL is DOUBLE PRECISION array, dimension (LDU, max(NN))
450 *> The (block lower triangular) left eigenvector matrix for
451 *> the matrices in S1 and P1. (See DTGEVC for the format.)
454 *> \param[out] EVECTR
456 *> EVECTR is DOUBLE PRECISION array, dimension (LDU, max(NN))
457 *> The (block upper triangular) right eigenvector matrix for
458 *> the matrices in S1 and P1. (See DTGEVC for the format.)
463 *> WORK is DOUBLE PRECISION array, dimension (LWORK)
469 *> The number of entries in WORK. This must be at least
470 *> max( 2 * N**2, 6*N, 1 ), for all N=NN(j).
473 *> \param[out] LLWORK
475 *> LLWORK is LOGICAL array, dimension (max(NN))
478 *> \param[out] RESULT
480 *> RESULT is DOUBLE PRECISION array, dimension (15)
481 *> The values computed by the tests described above.
482 *> The values are currently limited to 1/ulp, to avoid
489 *> = 0: successful exit
490 *> < 0: if INFO = -i, the i-th argument had an illegal value
491 *> > 0: A routine returned an error code. INFO is the
492 *> absolute value of the INFO value returned.
498 *> \author Univ. of Tennessee
499 *> \author Univ. of California Berkeley
500 *> \author Univ. of Colorado Denver
505 *> \ingroup double_eig
507 * =====================================================================
508 SUBROUTINE DCHKGG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
509 $ TSTDIF, THRSHN, NOUNIT, A, LDA, B, H, T, S1,
510 $ S2, P1, P2, U, LDU, V, Q, Z, ALPHR1, ALPHI1,
511 $ BETA1, ALPHR3, ALPHI3, BETA3, EVECTL, EVECTR,
512 $ WORK, LWORK, LLWORK, RESULT, INFO )
514 * -- LAPACK test routine (version 3.6.1) --
515 * -- LAPACK is a software package provided by Univ. of Tennessee, --
516 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
519 * .. Scalar Arguments ..
521 INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES
522 DOUBLE PRECISION THRESH, THRSHN
524 * .. Array Arguments ..
525 LOGICAL DOTYPE( * ), LLWORK( * )
526 INTEGER ISEED( 4 ), NN( * )
527 DOUBLE PRECISION A( LDA, * ), ALPHI1( * ), ALPHI3( * ),
528 $ ALPHR1( * ), ALPHR3( * ), B( LDA, * ),
529 $ BETA1( * ), BETA3( * ), EVECTL( LDU, * ),
530 $ EVECTR( LDU, * ), H( LDA, * ), P1( LDA, * ),
531 $ P2( LDA, * ), Q( LDU, * ), RESULT( 15 ),
532 $ S1( LDA, * ), S2( LDA, * ), T( LDA, * ),
533 $ U( LDU, * ), V( LDU, * ), WORK( * ),
537 * =====================================================================
540 DOUBLE PRECISION ZERO, ONE
541 PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0 )
543 PARAMETER ( MAXTYP = 26 )
545 * .. Local Scalars ..
547 INTEGER I1, IADD, IINFO, IN, J, JC, JR, JSIZE, JTYPE,
548 $ LWKOPT, MTYPES, N, N1, NERRS, NMATS, NMAX,
550 DOUBLE PRECISION ANORM, BNORM, SAFMAX, SAFMIN, TEMP1, TEMP2,
554 INTEGER IASIGN( MAXTYP ), IBSIGN( MAXTYP ),
555 $ IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ),
556 $ KATYPE( MAXTYP ), KAZERO( MAXTYP ),
557 $ KBMAGN( MAXTYP ), KBTYPE( MAXTYP ),
558 $ KBZERO( MAXTYP ), KCLASS( MAXTYP ),
559 $ KTRIAN( MAXTYP ), KZ1( 6 ), KZ2( 6 )
560 DOUBLE PRECISION DUMMA( 4 ), RMAGN( 0: 3 )
562 * .. External Functions ..
563 DOUBLE PRECISION DLAMCH, DLANGE, DLARND
564 EXTERNAL DLAMCH, DLANGE, DLARND
566 * .. External Subroutines ..
567 EXTERNAL DGEQR2, DGET51, DGET52, DGGHRD, DHGEQZ, DLABAD,
568 $ DLACPY, DLARFG, DLASET, DLASUM, DLATM4, DORM2R,
571 * .. Intrinsic Functions ..
572 INTRINSIC ABS, DBLE, MAX, MIN, SIGN
574 * .. Data statements ..
575 DATA KCLASS / 15*1, 10*2, 1*3 /
576 DATA KZ1 / 0, 1, 2, 1, 3, 3 /
577 DATA KZ2 / 0, 0, 1, 2, 1, 1 /
578 DATA KADD / 0, 0, 0, 0, 3, 2 /
579 DATA KATYPE / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
580 $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 /
581 DATA KBTYPE / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
582 $ 1, 1, -4, 2, -4, 8*8, 0 /
583 DATA KAZERO / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
585 DATA KBZERO / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
587 DATA KAMAGN / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
589 DATA KBMAGN / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
591 DATA KTRIAN / 16*0, 10*1 /
592 DATA IASIGN / 6*0, 2, 0, 2*2, 2*0, 3*2, 0, 2, 3*0,
594 DATA IBSIGN / 7*0, 2, 2*0, 2*2, 2*0, 2, 0, 2, 9*0 /
596 * .. Executable Statements ..
605 NMAX = MAX( NMAX, NN( J ) )
610 * Maximum blocksize and shift -- we assume that blocksize and number
611 * of shifts are monotone increasing functions of N.
613 LWKOPT = MAX( 6*NMAX, 2*NMAX*NMAX, 1 )
617 IF( NSIZES.LT.0 ) THEN
619 ELSE IF( BADNN ) THEN
621 ELSE IF( NTYPES.LT.0 ) THEN
623 ELSE IF( THRESH.LT.ZERO ) THEN
625 ELSE IF( LDA.LE.1 .OR. LDA.LT.NMAX ) THEN
627 ELSE IF( LDU.LE.1 .OR. LDU.LT.NMAX ) THEN
629 ELSE IF( LWKOPT.GT.LWORK ) THEN
634 CALL XERBLA( 'DCHKGG', -INFO )
638 * Quick return if possible
640 IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 )
643 SAFMIN = DLAMCH( 'Safe minimum' )
644 ULP = DLAMCH( 'Epsilon' )*DLAMCH( 'Base' )
645 SAFMIN = SAFMIN / ULP
646 SAFMAX = ONE / SAFMIN
647 CALL DLABAD( SAFMIN, SAFMAX )
650 * The values RMAGN(2:3) depend on N, see below.
655 * Loop over sizes, types
661 DO 240 JSIZE = 1, NSIZES
664 RMAGN( 2 ) = SAFMAX*ULP / DBLE( N1 )
665 RMAGN( 3 ) = SAFMIN*ULPINV*N1
667 IF( NSIZES.NE.1 ) THEN
668 MTYPES = MIN( MAXTYP, NTYPES )
670 MTYPES = MIN( MAXTYP+1, NTYPES )
673 DO 230 JTYPE = 1, MTYPES
674 IF( .NOT.DOTYPE( JTYPE ) )
679 * Save ISEED in case of an error.
682 IOLDSD( J ) = ISEED( J )
693 * Description of control parameters:
695 * KZLASS: =1 means w/o rotation, =2 means w/ rotation,
697 * KATYPE: the "type" to be passed to DLATM4 for computing A.
698 * KAZERO: the pattern of zeros on the diagonal for A:
699 * =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
700 * =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
701 * =6: ( 0, 1, 0, xxx, 0 ). (xxx means a string of
703 * KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
704 * =2: large, =3: small.
705 * IASIGN: 1 if the diagonal elements of A are to be
706 * multiplied by a random magnitude 1 number, =2 if
707 * randomly chosen diagonal blocks are to be rotated
708 * to form 2x2 blocks.
709 * KBTYPE, KBZERO, KBMAGN, IBSIGN: the same, but for B.
710 * KTRIAN: =0: don't fill in the upper triangle, =1: do.
711 * KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
712 * RMAGN: used to implement KAMAGN and KBMAGN.
714 IF( MTYPES.GT.MAXTYP )
717 IF( KCLASS( JTYPE ).LT.3 ) THEN
719 * Generate A (w/o rotation)
721 IF( ABS( KATYPE( JTYPE ) ).EQ.3 ) THEN
722 IN = 2*( ( N-1 ) / 2 ) + 1
724 $ CALL DLASET( 'Full', N, N, ZERO, ZERO, A, LDA )
728 CALL DLATM4( KATYPE( JTYPE ), IN, KZ1( KAZERO( JTYPE ) ),
729 $ KZ2( KAZERO( JTYPE ) ), IASIGN( JTYPE ),
730 $ RMAGN( KAMAGN( JTYPE ) ), ULP,
731 $ RMAGN( KTRIAN( JTYPE )*KAMAGN( JTYPE ) ), 2,
733 IADD = KADD( KAZERO( JTYPE ) )
734 IF( IADD.GT.0 .AND. IADD.LE.N )
735 $ A( IADD, IADD ) = RMAGN( KAMAGN( JTYPE ) )
737 * Generate B (w/o rotation)
739 IF( ABS( KBTYPE( JTYPE ) ).EQ.3 ) THEN
740 IN = 2*( ( N-1 ) / 2 ) + 1
742 $ CALL DLASET( 'Full', N, N, ZERO, ZERO, B, LDA )
746 CALL DLATM4( KBTYPE( JTYPE ), IN, KZ1( KBZERO( JTYPE ) ),
747 $ KZ2( KBZERO( JTYPE ) ), IBSIGN( JTYPE ),
748 $ RMAGN( KBMAGN( JTYPE ) ), ONE,
749 $ RMAGN( KTRIAN( JTYPE )*KBMAGN( JTYPE ) ), 2,
751 IADD = KADD( KBZERO( JTYPE ) )
752 IF( IADD.NE.0 .AND. IADD.LE.N )
753 $ B( IADD, IADD ) = RMAGN( KBMAGN( JTYPE ) )
755 IF( KCLASS( JTYPE ).EQ.2 .AND. N.GT.0 ) THEN
759 * Generate U, V as Householder transformations times
764 U( JR, JC ) = DLARND( 3, ISEED )
765 V( JR, JC ) = DLARND( 3, ISEED )
767 CALL DLARFG( N+1-JC, U( JC, JC ), U( JC+1, JC ), 1,
769 WORK( 2*N+JC ) = SIGN( ONE, U( JC, JC ) )
771 CALL DLARFG( N+1-JC, V( JC, JC ), V( JC+1, JC ), 1,
773 WORK( 3*N+JC ) = SIGN( ONE, V( JC, JC ) )
778 WORK( 3*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
781 WORK( 4*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
783 * Apply the diagonal matrices
787 A( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
789 B( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
793 CALL DORM2R( 'L', 'N', N, N, N-1, U, LDU, WORK, A,
794 $ LDA, WORK( 2*N+1 ), IINFO )
797 CALL DORM2R( 'R', 'T', N, N, N-1, V, LDU, WORK( N+1 ),
798 $ A, LDA, WORK( 2*N+1 ), IINFO )
801 CALL DORM2R( 'L', 'N', N, N, N-1, U, LDU, WORK, B,
802 $ LDA, WORK( 2*N+1 ), IINFO )
805 CALL DORM2R( 'R', 'T', N, N, N-1, V, LDU, WORK( N+1 ),
806 $ B, LDA, WORK( 2*N+1 ), IINFO )
816 A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
818 B( JR, JC ) = RMAGN( KBMAGN( JTYPE ) )*
824 ANORM = DLANGE( '1', N, N, A, LDA, WORK )
825 BNORM = DLANGE( '1', N, N, B, LDA, WORK )
829 IF( IINFO.NE.0 ) THEN
830 WRITE( NOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE,
838 * Call DGEQR2, DORM2R, and DGGHRD to compute H, T, U, and V
840 CALL DLACPY( ' ', N, N, A, LDA, H, LDA )
841 CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
845 CALL DGEQR2( N, N, T, LDA, WORK, WORK( N+1 ), IINFO )
846 IF( IINFO.NE.0 ) THEN
847 WRITE( NOUNIT, FMT = 9999 )'DGEQR2', IINFO, N, JTYPE,
853 CALL DORM2R( 'L', 'T', N, N, N, T, LDA, WORK, H, LDA,
854 $ WORK( N+1 ), IINFO )
855 IF( IINFO.NE.0 ) THEN
856 WRITE( NOUNIT, FMT = 9999 )'DORM2R', IINFO, N, JTYPE,
862 CALL DLASET( 'Full', N, N, ZERO, ONE, U, LDU )
863 CALL DORM2R( 'R', 'N', N, N, N, T, LDA, WORK, U, LDU,
864 $ WORK( N+1 ), IINFO )
865 IF( IINFO.NE.0 ) THEN
866 WRITE( NOUNIT, FMT = 9999 )'DORM2R', IINFO, N, JTYPE,
872 CALL DGGHRD( 'V', 'I', N, 1, N, H, LDA, T, LDA, U, LDU, V,
874 IF( IINFO.NE.0 ) THEN
875 WRITE( NOUNIT, FMT = 9999 )'DGGHRD', IINFO, N, JTYPE,
884 CALL DGET51( 1, N, A, LDA, H, LDA, U, LDU, V, LDU, WORK,
886 CALL DGET51( 1, N, B, LDA, T, LDA, U, LDU, V, LDU, WORK,
888 CALL DGET51( 3, N, B, LDA, T, LDA, U, LDU, U, LDU, WORK,
890 CALL DGET51( 3, N, B, LDA, T, LDA, V, LDU, V, LDU, WORK,
893 * Call DHGEQZ to compute S1, P1, S2, P2, Q, and Z, do tests.
899 CALL DLACPY( ' ', N, N, H, LDA, S2, LDA )
900 CALL DLACPY( ' ', N, N, T, LDA, P2, LDA )
904 CALL DHGEQZ( 'E', 'N', 'N', N, 1, N, S2, LDA, P2, LDA,
905 $ ALPHR3, ALPHI3, BETA3, Q, LDU, Z, LDU, WORK,
907 IF( IINFO.NE.0 ) THEN
908 WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(E)', IINFO, N, JTYPE,
914 * Eigenvalues and Full Schur Form
916 CALL DLACPY( ' ', N, N, H, LDA, S2, LDA )
917 CALL DLACPY( ' ', N, N, T, LDA, P2, LDA )
919 CALL DHGEQZ( 'S', 'N', 'N', N, 1, N, S2, LDA, P2, LDA,
920 $ ALPHR1, ALPHI1, BETA1, Q, LDU, Z, LDU, WORK,
922 IF( IINFO.NE.0 ) THEN
923 WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(S)', IINFO, N, JTYPE,
929 * Eigenvalues, Schur Form, and Schur Vectors
931 CALL DLACPY( ' ', N, N, H, LDA, S1, LDA )
932 CALL DLACPY( ' ', N, N, T, LDA, P1, LDA )
934 CALL DHGEQZ( 'S', 'I', 'I', N, 1, N, S1, LDA, P1, LDA,
935 $ ALPHR1, ALPHI1, BETA1, Q, LDU, Z, LDU, WORK,
937 IF( IINFO.NE.0 ) THEN
938 WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(V)', IINFO, N, JTYPE,
948 CALL DGET51( 1, N, H, LDA, S1, LDA, Q, LDU, Z, LDU, WORK,
950 CALL DGET51( 1, N, T, LDA, P1, LDA, Q, LDU, Z, LDU, WORK,
952 CALL DGET51( 3, N, T, LDA, P1, LDA, Q, LDU, Q, LDU, WORK,
954 CALL DGET51( 3, N, T, LDA, P1, LDA, Z, LDU, Z, LDU, WORK,
957 * Compute the Left and Right Eigenvectors of (S1,P1)
959 * 9: Compute the left eigenvector Matrix without
965 * To test "SELECT" option, compute half of the eigenvectors
966 * in one call, and half in another
973 LLWORK( J ) = .FALSE.
976 CALL DTGEVC( 'L', 'S', LLWORK, N, S1, LDA, P1, LDA, EVECTL,
977 $ LDU, DUMMA, LDU, N, IN, WORK, IINFO )
978 IF( IINFO.NE.0 ) THEN
979 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(L,S1)', IINFO, N,
987 LLWORK( J ) = .FALSE.
993 CALL DTGEVC( 'L', 'S', LLWORK, N, S1, LDA, P1, LDA,
994 $ EVECTL( 1, I1+1 ), LDU, DUMMA, LDU, N, IN,
996 IF( IINFO.NE.0 ) THEN
997 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(L,S2)', IINFO, N,
1003 CALL DGET52( .TRUE., N, S1, LDA, P1, LDA, EVECTL, LDU,
1004 $ ALPHR1, ALPHI1, BETA1, WORK, DUMMA( 1 ) )
1005 RESULT( 9 ) = DUMMA( 1 )
1006 IF( DUMMA( 2 ).GT.THRSHN ) THEN
1007 WRITE( NOUNIT, FMT = 9998 )'Left', 'DTGEVC(HOWMNY=S)',
1008 $ DUMMA( 2 ), N, JTYPE, IOLDSD
1011 * 10: Compute the left eigenvector Matrix with
1012 * back transforming:
1015 RESULT( 10 ) = ULPINV
1016 CALL DLACPY( 'F', N, N, Q, LDU, EVECTL, LDU )
1017 CALL DTGEVC( 'L', 'B', LLWORK, N, S1, LDA, P1, LDA, EVECTL,
1018 $ LDU, DUMMA, LDU, N, IN, WORK, IINFO )
1019 IF( IINFO.NE.0 ) THEN
1020 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(L,B)', IINFO, N,
1026 CALL DGET52( .TRUE., N, H, LDA, T, LDA, EVECTL, LDU, ALPHR1,
1027 $ ALPHI1, BETA1, WORK, DUMMA( 1 ) )
1028 RESULT( 10 ) = DUMMA( 1 )
1029 IF( DUMMA( 2 ).GT.THRSHN ) THEN
1030 WRITE( NOUNIT, FMT = 9998 )'Left', 'DTGEVC(HOWMNY=B)',
1031 $ DUMMA( 2 ), N, JTYPE, IOLDSD
1034 * 11: Compute the right eigenvector Matrix without
1035 * back transforming:
1038 RESULT( 11 ) = ULPINV
1040 * To test "SELECT" option, compute half of the eigenvectors
1041 * in one call, and half in another
1045 LLWORK( J ) = .TRUE.
1047 DO 170 J = I1 + 1, N
1048 LLWORK( J ) = .FALSE.
1051 CALL DTGEVC( 'R', 'S', LLWORK, N, S1, LDA, P1, LDA, DUMMA,
1052 $ LDU, EVECTR, LDU, N, IN, WORK, IINFO )
1053 IF( IINFO.NE.0 ) THEN
1054 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(R,S1)', IINFO, N,
1062 LLWORK( J ) = .FALSE.
1064 DO 190 J = I1 + 1, N
1065 LLWORK( J ) = .TRUE.
1068 CALL DTGEVC( 'R', 'S', LLWORK, N, S1, LDA, P1, LDA, DUMMA,
1069 $ LDU, EVECTR( 1, I1+1 ), LDU, N, IN, WORK,
1071 IF( IINFO.NE.0 ) THEN
1072 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(R,S2)', IINFO, N,
1078 CALL DGET52( .FALSE., N, S1, LDA, P1, LDA, EVECTR, LDU,
1079 $ ALPHR1, ALPHI1, BETA1, WORK, DUMMA( 1 ) )
1080 RESULT( 11 ) = DUMMA( 1 )
1081 IF( DUMMA( 2 ).GT.THRESH ) THEN
1082 WRITE( NOUNIT, FMT = 9998 )'Right', 'DTGEVC(HOWMNY=S)',
1083 $ DUMMA( 2 ), N, JTYPE, IOLDSD
1086 * 12: Compute the right eigenvector Matrix with
1087 * back transforming:
1090 RESULT( 12 ) = ULPINV
1091 CALL DLACPY( 'F', N, N, Z, LDU, EVECTR, LDU )
1092 CALL DTGEVC( 'R', 'B', LLWORK, N, S1, LDA, P1, LDA, DUMMA,
1093 $ LDU, EVECTR, LDU, N, IN, WORK, IINFO )
1094 IF( IINFO.NE.0 ) THEN
1095 WRITE( NOUNIT, FMT = 9999 )'DTGEVC(R,B)', IINFO, N,
1101 CALL DGET52( .FALSE., N, H, LDA, T, LDA, EVECTR, LDU,
1102 $ ALPHR1, ALPHI1, BETA1, WORK, DUMMA( 1 ) )
1103 RESULT( 12 ) = DUMMA( 1 )
1104 IF( DUMMA( 2 ).GT.THRESH ) THEN
1105 WRITE( NOUNIT, FMT = 9998 )'Right', 'DTGEVC(HOWMNY=B)',
1106 $ DUMMA( 2 ), N, JTYPE, IOLDSD
1109 * Tests 13--15 are done only on request
1115 CALL DGET51( 2, N, S1, LDA, S2, LDA, Q, LDU, Z, LDU,
1116 $ WORK, RESULT( 13 ) )
1117 CALL DGET51( 2, N, P1, LDA, P2, LDA, Q, LDU, Z, LDU,
1118 $ WORK, RESULT( 14 ) )
1125 TEMP1 = MAX( TEMP1, ABS( ALPHR1( J )-ALPHR3( J ) )+
1126 $ ABS( ALPHI1( J )-ALPHI3( J ) ) )
1127 TEMP2 = MAX( TEMP2, ABS( BETA1( J )-BETA3( J ) ) )
1130 TEMP1 = TEMP1 / MAX( SAFMIN, ULP*MAX( TEMP1, ANORM ) )
1131 TEMP2 = TEMP2 / MAX( SAFMIN, ULP*MAX( TEMP2, BNORM ) )
1132 RESULT( 15 ) = MAX( TEMP1, TEMP2 )
1141 * End of Loop -- Check for RESULT(j) > THRESH
1145 NTESTT = NTESTT + NTEST
1147 * Print out tests which fail.
1149 DO 220 JR = 1, NTEST
1150 IF( RESULT( JR ).GE.THRESH ) THEN
1152 * If this is the first test to fail,
1153 * print a header to the data file.
1155 IF( NERRS.EQ.0 ) THEN
1156 WRITE( NOUNIT, FMT = 9997 )'DGG'
1160 WRITE( NOUNIT, FMT = 9996 )
1161 WRITE( NOUNIT, FMT = 9995 )
1162 WRITE( NOUNIT, FMT = 9994 )'Orthogonal'
1166 WRITE( NOUNIT, FMT = 9993 )'orthogonal', '''',
1167 $ 'transpose', ( '''', J = 1, 10 )
1171 IF( RESULT( JR ).LT.10000.0D0 ) THEN
1172 WRITE( NOUNIT, FMT = 9992 )N, JTYPE, IOLDSD, JR,
1175 WRITE( NOUNIT, FMT = 9991 )N, JTYPE, IOLDSD, JR,
1186 CALL DLASUM( 'DGG', NOUNIT, NERRS, NTESTT )
1189 9999 FORMAT( ' DCHKGG: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
1190 $ I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' )
1192 9998 FORMAT( ' DCHKGG: ', A, ' Eigenvectors from ', A, ' incorrectly ',
1193 $ 'normalized.', / ' Bits of error=', 0P, G10.3, ',', 9X,
1194 $ 'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5,
1197 9997 FORMAT( / 1X, A3, ' -- Real Generalized eigenvalue problem' )
1199 9996 FORMAT( ' Matrix types (see DCHKGG for details): ' )
1201 9995 FORMAT( ' Special Matrices:', 23X,
1202 $ '(J''=transposed Jordan block)',
1203 $ / ' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ',
1204 $ '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices: ( ',
1205 $ 'D=diag(0,1,2,...) )', / ' 7=(D,I) 9=(large*D, small*I',
1206 $ ') 11=(large*I, small*D) 13=(large*D, large*I)', /
1207 $ ' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ',
1208 $ ' 14=(small*D, small*I)', / ' 15=(D, reversed D)' )
1209 9994 FORMAT( ' Matrices Rotated by Random ', A, ' Matrices U, V:',
1210 $ / ' 16=Transposed Jordan Blocks 19=geometric ',
1211 $ 'alpha, beta=0,1', / ' 17=arithm. alpha&beta ',
1212 $ ' 20=arithmetic alpha, beta=0,1', / ' 18=clustered ',
1213 $ 'alpha, beta=0,1 21=random alpha, beta=0,1',
1214 $ / ' Large & Small Matrices:', / ' 22=(large, small) ',
1215 $ '23=(small,large) 24=(small,small) 25=(large,large)',
1216 $ / ' 26=random O(1) matrices.' )
1218 9993 FORMAT( / ' Tests performed: (H is Hessenberg, S is Schur, B, ',
1219 $ 'T, P are triangular,', / 20X, 'U, V, Q, and Z are ', A,
1220 $ ', l and r are the', / 20X,
1221 $ 'appropriate left and right eigenvectors, resp., a is',
1222 $ / 20X, 'alpha, b is beta, and ', A, ' means ', A, '.)',
1223 $ / ' 1 = | A - U H V', A,
1224 $ ' | / ( |A| n ulp ) 2 = | B - U T V', A,
1225 $ ' | / ( |B| n ulp )', / ' 3 = | I - UU', A,
1226 $ ' | / ( n ulp ) 4 = | I - VV', A,
1227 $ ' | / ( n ulp )', / ' 5 = | H - Q S Z', A,
1228 $ ' | / ( |H| n ulp )', 6X, '6 = | T - Q P Z', A,
1229 $ ' | / ( |T| n ulp )', / ' 7 = | I - QQ', A,
1230 $ ' | / ( n ulp ) 8 = | I - ZZ', A,
1231 $ ' | / ( n ulp )', / ' 9 = max | ( b S - a P )', A,
1232 $ ' l | / const. 10 = max | ( b H - a T )', A,
1233 $ ' l | / const.', /
1234 $ ' 11= max | ( b S - a P ) r | / const. 12 = max | ( b H',
1235 $ ' - a T ) r | / const.', / 1X )
1237 9992 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
1238 $ 4( I4, ',' ), ' result ', I2, ' is', 0P, F8.2 )
1239 9991 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
1240 $ 4( I4, ',' ), ' result ', I2, ' is', 1P, D10.3 )