STYLE: Remove trailing whitespace in Fortran files
[platform/upstream/lapack.git] / TESTING / EIG / dchkgg.f
1 *> \brief \b DCHKGG
2 *
3 *  =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 *            http://www.netlib.org/lapack/explore-html/
7 *
8 *  Definition:
9 *  ===========
10 *
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 )
16 *
17 *       .. Scalar Arguments ..
18 *       LOGICAL            TSTDIF
19 *       INTEGER            INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES
20 *       DOUBLE PRECISION   THRESH, THRSHN
21 *       ..
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( * ),
32 *      $                   Z( LDU, * )
33 *       ..
34 *
35 *
36 *> \par Purpose:
37 *  =============
38 *>
39 *> \verbatim
40 *>
41 *> DCHKGG  checks the nonsymmetric generalized eigenvalue problem
42 *> routines.
43 *>                                T          T        T
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
46 *> orthogonal.
47 *>                                 T          T
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
56 *> problem
57 *>
58 *>     det( A - w(j) B ) = 0
59 *>
60 *> and m(j) = beta(j)/alpha(j) is a root of the essentially equivalent
61 *> problem
62 *>
63 *>     det( m(j) A - B ) = 0
64 *>
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).
69 *>
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:
76 *>
77 *>                  T
78 *> (1)   | A - U H V  | / ( |A| n ulp )
79 *>
80 *>                  T
81 *> (2)   | B - U T V  | / ( |B| n ulp )
82 *>
83 *>               T
84 *> (3)   | I - UU  | / ( n ulp )
85 *>
86 *>               T
87 *> (4)   | I - VV  | / ( n ulp )
88 *>
89 *>                  T
90 *> (5)   | H - Q S Z  | / ( |H| n ulp )
91 *>
92 *>                  T
93 *> (6)   | T - Q P Z  | / ( |T| n ulp )
94 *>
95 *>               T
96 *> (7)   | I - QQ  | / ( n ulp )
97 *>
98 *>               T
99 *> (8)   | I - ZZ  | / ( n ulp )
100 *>
101 *> (9)   max over all left eigenvalue/-vector pairs (beta/alpha,l) of
102 *>
103 *>    | l**H * (beta S - alpha P) | / ( ulp max( |beta S|, |alpha P| ) )
104 *>
105 *> (10)  max over all left eigenvalue/-vector pairs (beta/alpha,l') of
106 *>                           T
107 *>   | l'**H * (beta H - alpha T) | / ( ulp max( |beta H|, |alpha T| ) )
108 *>
109 *>       where the eigenvectors l' are the result of passing Q to
110 *>       DTGEVC and back transforming (HOWMNY='B').
111 *>
112 *> (11)  max over all right eigenvalue/-vector pairs (beta/alpha,r) of
113 *>
114 *>       | (beta S - alpha T) r | / ( ulp max( |beta S|, |alpha T| ) )
115 *>
116 *> (12)  max over all right eigenvalue/-vector pairs (beta/alpha,r') of
117 *>
118 *>       | (beta H - alpha T) r' | / ( ulp max( |beta H|, |alpha T| ) )
119 *>
120 *>       where the eigenvectors r' are the result of passing Z to
121 *>       DTGEVC and back transforming (HOWMNY='B').
122 *>
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.
126 *>
127 *> (13)  | S(Q,Z computed) - S(Q,Z not computed) | / ( |S| ulp )
128 *>
129 *> (14)  | P(Q,Z computed) - P(Q,Z not computed) | / ( |P| ulp )
130 *>
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
133 *>
134 *> In addition, the normalization of L and R are checked, and compared
135 *> with the threshold THRSHN.
136 *>
137 *> Test Matrices
138 *> ---- --------
139 *>
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:
145 *>
146 *> (1)  ( 0, 0 )         (a pair of zero matrices)
147 *>
148 *> (2)  ( I, 0 )         (an identity and a zero matrix)
149 *>
150 *> (3)  ( 0, I )         (an identity and a zero matrix)
151 *>
152 *> (4)  ( I, I )         (a pair of identity matrices)
153 *>
154 *>         t   t
155 *> (5)  ( J , J  )       (a pair of transposed Jordan blocks)
156 *>
157 *>                                     t                ( I   0  )
158 *> (6)  ( X, Y )         where  X = ( J   0  )  and Y = (      t )
159 *>                                  ( 0   I  )          ( 0   J  )
160 *>                       and I is a k x k identity and J a (k+1)x(k+1)
161 *>                       Jordan block; k=(N-1)/2
162 *>
163 *> (7)  ( D, I )         where D is diag( 0, 1,..., N-1 ) (a diagonal
164 *>                       matrix with those diagonal entries.)
165 *> (8)  ( I, D )
166 *>
167 *> (9)  ( big*D, small*I ) where "big" is near overflow and small=1/big
168 *>
169 *> (10) ( small*D, big*I )
170 *>
171 *> (11) ( big*I, small*D )
172 *>
173 *> (12) ( small*I, big*D )
174 *>
175 *> (13) ( big*D, big*I )
176 *>
177 *> (14) ( small*D, small*I )
178 *>
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 )
181 *>           t   t
182 *> (16) U ( J , J ) V     where U and V are random orthogonal matrices.
183 *>
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 )
189 *>
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.
193 *>
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 )
196 *>
197 *>                                                        N-5
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 )
200 *>
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.
204 *>
205 *> (22) U ( big*T1, small*T2 ) V    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
206 *>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
207 *>
208 *> (23) U ( small*T1, big*T2 ) V    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
209 *>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
210 *>
211 *> (24) U ( small*T1, small*T2 ) V  diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
212 *>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
213 *>
214 *> (25) U ( big*T1, big*T2 ) V      diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
215 *>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
216 *>
217 *> (26) U ( T1, T2 ) V     where T1 and T2 are random upper-triangular
218 *>                         matrices.
219 *> \endverbatim
220 *
221 *  Arguments:
222 *  ==========
223 *
224 *> \param[in] NSIZES
225 *> \verbatim
226 *>          NSIZES is INTEGER
227 *>          The number of sizes of matrices to use.  If it is zero,
228 *>          DCHKGG does nothing.  It must be at least zero.
229 *> \endverbatim
230 *>
231 *> \param[in] NN
232 *> \verbatim
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
236 *>          zero.
237 *> \endverbatim
238 *>
239 *> \param[in] NTYPES
240 *> \verbatim
241 *>          NTYPES is INTEGER
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. .
248 *> \endverbatim
249 *>
250 *> \param[in] DOTYPE
251 *> \verbatim
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)
259 *>          will be ignored.
260 *> \endverbatim
261 *>
262 *> \param[in,out] ISEED
263 *> \verbatim
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
273 *>          sequence.
274 *> \endverbatim
275 *>
276 *> \param[in] THRESH
277 *> \verbatim
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.
285 *> \endverbatim
286 *>
287 *> \param[in] TSTDIF
288 *> \verbatim
289 *>          TSTDIF is LOGICAL
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
295 *>                     tested.
296 *> \endverbatim
297 *>
298 *> \param[in] THRSHN
299 *> \verbatim
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.)
307 *> \endverbatim
308 *>
309 *> \param[in] NOUNIT
310 *> \verbatim
311 *>          NOUNIT is INTEGER
312 *>          The FORTRAN unit number for printing out error messages
313 *>          (e.g., if a routine returns IINFO not equal to 0.)
314 *> \endverbatim
315 *>
316 *> \param[in,out] A
317 *> \verbatim
318 *>          A is DOUBLE PRECISION array, dimension
319 *>                            (LDA, max(NN))
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.
323 *> \endverbatim
324 *>
325 *> \param[in] LDA
326 *> \verbatim
327 *>          LDA is INTEGER
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 ).
330 *> \endverbatim
331 *>
332 *> \param[in,out] B
333 *> \verbatim
334 *>          B is DOUBLE PRECISION array, dimension
335 *>                            (LDA, max(NN))
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.
339 *> \endverbatim
340 *>
341 *> \param[out] H
342 *> \verbatim
343 *>          H is DOUBLE PRECISION array, dimension (LDA, max(NN))
344 *>          The upper Hessenberg matrix computed from A by DGGHRD.
345 *> \endverbatim
346 *>
347 *> \param[out] T
348 *> \verbatim
349 *>          T is DOUBLE PRECISION array, dimension (LDA, max(NN))
350 *>          The upper triangular matrix computed from B by DGGHRD.
351 *> \endverbatim
352 *>
353 *> \param[out] S1
354 *> \verbatim
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.
358 *> \endverbatim
359 *>
360 *> \param[out] S2
361 *> \verbatim
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.
365 *> \endverbatim
366 *>
367 *> \param[out] P1
368 *> \verbatim
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.
372 *> \endverbatim
373 *>
374 *> \param[out] P2
375 *> \verbatim
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.
379 *> \endverbatim
380 *>
381 *> \param[out] U
382 *> \verbatim
383 *>          U is DOUBLE PRECISION array, dimension (LDU, max(NN))
384 *>          The (left) orthogonal matrix computed by DGGHRD.
385 *> \endverbatim
386 *>
387 *> \param[in] LDU
388 *> \verbatim
389 *>          LDU is INTEGER
390 *>          The leading dimension of U, V, Q, Z, EVECTL, and EVEZTR.  It
391 *>          must be at least 1 and at least max( NN ).
392 *> \endverbatim
393 *>
394 *> \param[out] V
395 *> \verbatim
396 *>          V is DOUBLE PRECISION array, dimension (LDU, max(NN))
397 *>          The (right) orthogonal matrix computed by DGGHRD.
398 *> \endverbatim
399 *>
400 *> \param[out] Q
401 *> \verbatim
402 *>          Q is DOUBLE PRECISION array, dimension (LDU, max(NN))
403 *>          The (left) orthogonal matrix computed by DHGEQZ.
404 *> \endverbatim
405 *>
406 *> \param[out] Z
407 *> \verbatim
408 *>          Z is DOUBLE PRECISION array, dimension (LDU, max(NN))
409 *>          The (left) orthogonal matrix computed by DHGEQZ.
410 *> \endverbatim
411 *>
412 *> \param[out] ALPHR1
413 *> \verbatim
414 *>          ALPHR1 is DOUBLE PRECISION array, dimension (max(NN))
415 *> \endverbatim
416 *>
417 *> \param[out] ALPHI1
418 *> \verbatim
419 *>          ALPHI1 is DOUBLE PRECISION array, dimension (max(NN))
420 *> \endverbatim
421 *>
422 *> \param[out] BETA1
423 *> \verbatim
424 *>          BETA1 is DOUBLE PRECISION array, dimension (max(NN))
425 *>
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.
430 *> \endverbatim
431 *>
432 *> \param[out] ALPHR3
433 *> \verbatim
434 *>          ALPHR3 is DOUBLE PRECISION array, dimension (max(NN))
435 *> \endverbatim
436 *>
437 *> \param[out] ALPHI3
438 *> \verbatim
439 *>          ALPHI3 is DOUBLE PRECISION array, dimension (max(NN))
440 *> \endverbatim
441 *>
442 *> \param[out] BETA3
443 *> \verbatim
444 *>          BETA3 is DOUBLE PRECISION array, dimension (max(NN))
445 *> \endverbatim
446 *>
447 *> \param[out] EVECTL
448 *> \verbatim
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.)
452 *> \endverbatim
453 *>
454 *> \param[out] EVECTR
455 *> \verbatim
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.)
459 *> \endverbatim
460 *>
461 *> \param[out] WORK
462 *> \verbatim
463 *>          WORK is DOUBLE PRECISION array, dimension (LWORK)
464 *> \endverbatim
465 *>
466 *> \param[in] LWORK
467 *> \verbatim
468 *>          LWORK is INTEGER
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).
471 *> \endverbatim
472 *>
473 *> \param[out] LLWORK
474 *> \verbatim
475 *>          LLWORK is LOGICAL array, dimension (max(NN))
476 *> \endverbatim
477 *>
478 *> \param[out] RESULT
479 *> \verbatim
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
483 *>          overflow.
484 *> \endverbatim
485 *>
486 *> \param[out] INFO
487 *> \verbatim
488 *>          INFO is INTEGER
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.
493 *> \endverbatim
494 *
495 *  Authors:
496 *  ========
497 *
498 *> \author Univ. of Tennessee
499 *> \author Univ. of California Berkeley
500 *> \author Univ. of Colorado Denver
501 *> \author NAG Ltd.
502 *
503 *> \date June 2016
504 *
505 *> \ingroup double_eig
506 *
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 )
513 *
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..--
517 *     June 2016
518 *
519 *     .. Scalar Arguments ..
520       LOGICAL            TSTDIF
521       INTEGER            INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES
522       DOUBLE PRECISION   THRESH, THRSHN
523 *     ..
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( * ),
534      $                   Z( LDU, * )
535 *     ..
536 *
537 *  =====================================================================
538 *
539 *     .. Parameters ..
540       DOUBLE PRECISION   ZERO, ONE
541       PARAMETER          ( ZERO = 0.0D0, ONE = 1.0D0 )
542       INTEGER            MAXTYP
543       PARAMETER          ( MAXTYP = 26 )
544 *     ..
545 *     .. Local Scalars ..
546       LOGICAL            BADNN
547       INTEGER            I1, IADD, IINFO, IN, J, JC, JR, JSIZE, JTYPE,
548      $                   LWKOPT, MTYPES, N, N1, NERRS, NMATS, NMAX,
549      $                   NTEST, NTESTT
550       DOUBLE PRECISION   ANORM, BNORM, SAFMAX, SAFMIN, TEMP1, TEMP2,
551      $                   ULP, ULPINV
552 *     ..
553 *     .. Local Arrays ..
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 )
561 *     ..
562 *     .. External Functions ..
563       DOUBLE PRECISION   DLAMCH, DLANGE, DLARND
564       EXTERNAL           DLAMCH, DLANGE, DLARND
565 *     ..
566 *     .. External Subroutines ..
567       EXTERNAL           DGEQR2, DGET51, DGET52, DGGHRD, DHGEQZ, DLABAD,
568      $                   DLACPY, DLARFG, DLASET, DLASUM, DLATM4, DORM2R,
569      $                   DTGEVC, XERBLA
570 *     ..
571 *     .. Intrinsic Functions ..
572       INTRINSIC          ABS, DBLE, MAX, MIN, SIGN
573 *     ..
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,
584      $                   4*5, 4*3, 1 /
585       DATA               KBZERO / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
586      $                   4*6, 4*4, 1 /
587       DATA               KAMAGN / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
588      $                   2, 1 /
589       DATA               KBMAGN / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
590      $                   2, 1 /
591       DATA               KTRIAN / 16*0, 10*1 /
592       DATA               IASIGN / 6*0, 2, 0, 2*2, 2*0, 3*2, 0, 2, 3*0,
593      $                   5*2, 0 /
594       DATA               IBSIGN / 7*0, 2, 2*0, 2*2, 2*0, 2, 0, 2, 9*0 /
595 *     ..
596 *     .. Executable Statements ..
597 *
598 *     Check for errors
599 *
600       INFO = 0
601 *
602       BADNN = .FALSE.
603       NMAX = 1
604       DO 10 J = 1, NSIZES
605          NMAX = MAX( NMAX, NN( J ) )
606          IF( NN( J ).LT.0 )
607      $      BADNN = .TRUE.
608    10 CONTINUE
609 *
610 *     Maximum blocksize and shift -- we assume that blocksize and number
611 *     of shifts are monotone increasing functions of N.
612 *
613       LWKOPT = MAX( 6*NMAX, 2*NMAX*NMAX, 1 )
614 *
615 *     Check for errors
616 *
617       IF( NSIZES.LT.0 ) THEN
618          INFO = -1
619       ELSE IF( BADNN ) THEN
620          INFO = -2
621       ELSE IF( NTYPES.LT.0 ) THEN
622          INFO = -3
623       ELSE IF( THRESH.LT.ZERO ) THEN
624          INFO = -6
625       ELSE IF( LDA.LE.1 .OR. LDA.LT.NMAX ) THEN
626          INFO = -10
627       ELSE IF( LDU.LE.1 .OR. LDU.LT.NMAX ) THEN
628          INFO = -19
629       ELSE IF( LWKOPT.GT.LWORK ) THEN
630          INFO = -30
631       END IF
632 *
633       IF( INFO.NE.0 ) THEN
634          CALL XERBLA( 'DCHKGG', -INFO )
635          RETURN
636       END IF
637 *
638 *     Quick return if possible
639 *
640       IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 )
641      $   RETURN
642 *
643       SAFMIN = DLAMCH( 'Safe minimum' )
644       ULP = DLAMCH( 'Epsilon' )*DLAMCH( 'Base' )
645       SAFMIN = SAFMIN / ULP
646       SAFMAX = ONE / SAFMIN
647       CALL DLABAD( SAFMIN, SAFMAX )
648       ULPINV = ONE / ULP
649 *
650 *     The values RMAGN(2:3) depend on N, see below.
651 *
652       RMAGN( 0 ) = ZERO
653       RMAGN( 1 ) = ONE
654 *
655 *     Loop over sizes, types
656 *
657       NTESTT = 0
658       NERRS = 0
659       NMATS = 0
660 *
661       DO 240 JSIZE = 1, NSIZES
662          N = NN( JSIZE )
663          N1 = MAX( 1, N )
664          RMAGN( 2 ) = SAFMAX*ULP / DBLE( N1 )
665          RMAGN( 3 ) = SAFMIN*ULPINV*N1
666 *
667          IF( NSIZES.NE.1 ) THEN
668             MTYPES = MIN( MAXTYP, NTYPES )
669          ELSE
670             MTYPES = MIN( MAXTYP+1, NTYPES )
671          END IF
672 *
673          DO 230 JTYPE = 1, MTYPES
674             IF( .NOT.DOTYPE( JTYPE ) )
675      $         GO TO 230
676             NMATS = NMATS + 1
677             NTEST = 0
678 *
679 *           Save ISEED in case of an error.
680 *
681             DO 20 J = 1, 4
682                IOLDSD( J ) = ISEED( J )
683    20       CONTINUE
684 *
685 *           Initialize RESULT
686 *
687             DO 30 J = 1, 15
688                RESULT( J ) = ZERO
689    30       CONTINUE
690 *
691 *           Compute A and B
692 *
693 *           Description of control parameters:
694 *
695 *           KZLASS: =1 means w/o rotation, =2 means w/ rotation,
696 *                   =3 means random.
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
702 *                   non-zero entries.)
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.
713 *
714             IF( MTYPES.GT.MAXTYP )
715      $         GO TO 110
716             IINFO = 0
717             IF( KCLASS( JTYPE ).LT.3 ) THEN
718 *
719 *              Generate A (w/o rotation)
720 *
721                IF( ABS( KATYPE( JTYPE ) ).EQ.3 ) THEN
722                   IN = 2*( ( N-1 ) / 2 ) + 1
723                   IF( IN.NE.N )
724      $               CALL DLASET( 'Full', N, N, ZERO, ZERO, A, LDA )
725                ELSE
726                   IN = N
727                END IF
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,
732      $                      ISEED, A, LDA )
733                IADD = KADD( KAZERO( JTYPE ) )
734                IF( IADD.GT.0 .AND. IADD.LE.N )
735      $            A( IADD, IADD ) = RMAGN( KAMAGN( JTYPE ) )
736 *
737 *              Generate B (w/o rotation)
738 *
739                IF( ABS( KBTYPE( JTYPE ) ).EQ.3 ) THEN
740                   IN = 2*( ( N-1 ) / 2 ) + 1
741                   IF( IN.NE.N )
742      $               CALL DLASET( 'Full', N, N, ZERO, ZERO, B, LDA )
743                ELSE
744                   IN = N
745                END IF
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,
750      $                      ISEED, B, LDA )
751                IADD = KADD( KBZERO( JTYPE ) )
752                IF( IADD.NE.0 .AND. IADD.LE.N )
753      $            B( IADD, IADD ) = RMAGN( KBMAGN( JTYPE ) )
754 *
755                IF( KCLASS( JTYPE ).EQ.2 .AND. N.GT.0 ) THEN
756 *
757 *                 Include rotations
758 *
759 *                 Generate U, V as Householder transformations times
760 *                 a diagonal matrix.
761 *
762                   DO 50 JC = 1, N - 1
763                      DO 40 JR = JC, N
764                         U( JR, JC ) = DLARND( 3, ISEED )
765                         V( JR, JC ) = DLARND( 3, ISEED )
766    40                CONTINUE
767                      CALL DLARFG( N+1-JC, U( JC, JC ), U( JC+1, JC ), 1,
768      $                            WORK( JC ) )
769                      WORK( 2*N+JC ) = SIGN( ONE, U( JC, JC ) )
770                      U( JC, JC ) = ONE
771                      CALL DLARFG( N+1-JC, V( JC, JC ), V( JC+1, JC ), 1,
772      $                            WORK( N+JC ) )
773                      WORK( 3*N+JC ) = SIGN( ONE, V( JC, JC ) )
774                      V( JC, JC ) = ONE
775    50             CONTINUE
776                   U( N, N ) = ONE
777                   WORK( N ) = ZERO
778                   WORK( 3*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
779                   V( N, N ) = ONE
780                   WORK( 2*N ) = ZERO
781                   WORK( 4*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
782 *
783 *                 Apply the diagonal matrices
784 *
785                   DO 70 JC = 1, N
786                      DO 60 JR = 1, N
787                         A( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
788      $                                A( JR, JC )
789                         B( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
790      $                                B( JR, JC )
791    60                CONTINUE
792    70             CONTINUE
793                   CALL DORM2R( 'L', 'N', N, N, N-1, U, LDU, WORK, A,
794      $                         LDA, WORK( 2*N+1 ), IINFO )
795                   IF( IINFO.NE.0 )
796      $               GO TO 100
797                   CALL DORM2R( 'R', 'T', N, N, N-1, V, LDU, WORK( N+1 ),
798      $                         A, LDA, WORK( 2*N+1 ), IINFO )
799                   IF( IINFO.NE.0 )
800      $               GO TO 100
801                   CALL DORM2R( 'L', 'N', N, N, N-1, U, LDU, WORK, B,
802      $                         LDA, WORK( 2*N+1 ), IINFO )
803                   IF( IINFO.NE.0 )
804      $               GO TO 100
805                   CALL DORM2R( 'R', 'T', N, N, N-1, V, LDU, WORK( N+1 ),
806      $                         B, LDA, WORK( 2*N+1 ), IINFO )
807                   IF( IINFO.NE.0 )
808      $               GO TO 100
809                END IF
810             ELSE
811 *
812 *              Random matrices
813 *
814                DO 90 JC = 1, N
815                   DO 80 JR = 1, N
816                      A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
817      $                             DLARND( 2, ISEED )
818                      B( JR, JC ) = RMAGN( KBMAGN( JTYPE ) )*
819      $                             DLARND( 2, ISEED )
820    80             CONTINUE
821    90          CONTINUE
822             END IF
823 *
824             ANORM = DLANGE( '1', N, N, A, LDA, WORK )
825             BNORM = DLANGE( '1', N, N, B, LDA, WORK )
826 *
827   100       CONTINUE
828 *
829             IF( IINFO.NE.0 ) THEN
830                WRITE( NOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE,
831      $            IOLDSD
832                INFO = ABS( IINFO )
833                RETURN
834             END IF
835 *
836   110       CONTINUE
837 *
838 *           Call DGEQR2, DORM2R, and DGGHRD to compute H, T, U, and V
839 *
840             CALL DLACPY( ' ', N, N, A, LDA, H, LDA )
841             CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
842             NTEST = 1
843             RESULT( 1 ) = ULPINV
844 *
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,
848      $            IOLDSD
849                INFO = ABS( IINFO )
850                GO TO 210
851             END IF
852 *
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,
857      $            IOLDSD
858                INFO = ABS( IINFO )
859                GO TO 210
860             END IF
861 *
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,
867      $            IOLDSD
868                INFO = ABS( IINFO )
869                GO TO 210
870             END IF
871 *
872             CALL DGGHRD( 'V', 'I', N, 1, N, H, LDA, T, LDA, U, LDU, V,
873      $                   LDU, IINFO )
874             IF( IINFO.NE.0 ) THEN
875                WRITE( NOUNIT, FMT = 9999 )'DGGHRD', IINFO, N, JTYPE,
876      $            IOLDSD
877                INFO = ABS( IINFO )
878                GO TO 210
879             END IF
880             NTEST = 4
881 *
882 *           Do tests 1--4
883 *
884             CALL DGET51( 1, N, A, LDA, H, LDA, U, LDU, V, LDU, WORK,
885      $                   RESULT( 1 ) )
886             CALL DGET51( 1, N, B, LDA, T, LDA, U, LDU, V, LDU, WORK,
887      $                   RESULT( 2 ) )
888             CALL DGET51( 3, N, B, LDA, T, LDA, U, LDU, U, LDU, WORK,
889      $                   RESULT( 3 ) )
890             CALL DGET51( 3, N, B, LDA, T, LDA, V, LDU, V, LDU, WORK,
891      $                   RESULT( 4 ) )
892 *
893 *           Call DHGEQZ to compute S1, P1, S2, P2, Q, and Z, do tests.
894 *
895 *           Compute T1 and UZ
896 *
897 *           Eigenvalues only
898 *
899             CALL DLACPY( ' ', N, N, H, LDA, S2, LDA )
900             CALL DLACPY( ' ', N, N, T, LDA, P2, LDA )
901             NTEST = 5
902             RESULT( 5 ) = ULPINV
903 *
904             CALL DHGEQZ( 'E', 'N', 'N', N, 1, N, S2, LDA, P2, LDA,
905      $                   ALPHR3, ALPHI3, BETA3, Q, LDU, Z, LDU, WORK,
906      $                   LWORK, IINFO )
907             IF( IINFO.NE.0 ) THEN
908                WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(E)', IINFO, N, JTYPE,
909      $            IOLDSD
910                INFO = ABS( IINFO )
911                GO TO 210
912             END IF
913 *
914 *           Eigenvalues and Full Schur Form
915 *
916             CALL DLACPY( ' ', N, N, H, LDA, S2, LDA )
917             CALL DLACPY( ' ', N, N, T, LDA, P2, LDA )
918 *
919             CALL DHGEQZ( 'S', 'N', 'N', N, 1, N, S2, LDA, P2, LDA,
920      $                   ALPHR1, ALPHI1, BETA1, Q, LDU, Z, LDU, WORK,
921      $                   LWORK, IINFO )
922             IF( IINFO.NE.0 ) THEN
923                WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(S)', IINFO, N, JTYPE,
924      $            IOLDSD
925                INFO = ABS( IINFO )
926                GO TO 210
927             END IF
928 *
929 *           Eigenvalues, Schur Form, and Schur Vectors
930 *
931             CALL DLACPY( ' ', N, N, H, LDA, S1, LDA )
932             CALL DLACPY( ' ', N, N, T, LDA, P1, LDA )
933 *
934             CALL DHGEQZ( 'S', 'I', 'I', N, 1, N, S1, LDA, P1, LDA,
935      $                   ALPHR1, ALPHI1, BETA1, Q, LDU, Z, LDU, WORK,
936      $                   LWORK, IINFO )
937             IF( IINFO.NE.0 ) THEN
938                WRITE( NOUNIT, FMT = 9999 )'DHGEQZ(V)', IINFO, N, JTYPE,
939      $            IOLDSD
940                INFO = ABS( IINFO )
941                GO TO 210
942             END IF
943 *
944             NTEST = 8
945 *
946 *           Do Tests 5--8
947 *
948             CALL DGET51( 1, N, H, LDA, S1, LDA, Q, LDU, Z, LDU, WORK,
949      $                   RESULT( 5 ) )
950             CALL DGET51( 1, N, T, LDA, P1, LDA, Q, LDU, Z, LDU, WORK,
951      $                   RESULT( 6 ) )
952             CALL DGET51( 3, N, T, LDA, P1, LDA, Q, LDU, Q, LDU, WORK,
953      $                   RESULT( 7 ) )
954             CALL DGET51( 3, N, T, LDA, P1, LDA, Z, LDU, Z, LDU, WORK,
955      $                   RESULT( 8 ) )
956 *
957 *           Compute the Left and Right Eigenvectors of (S1,P1)
958 *
959 *           9: Compute the left eigenvector Matrix without
960 *              back transforming:
961 *
962             NTEST = 9
963             RESULT( 9 ) = ULPINV
964 *
965 *           To test "SELECT" option, compute half of the eigenvectors
966 *           in one call, and half in another
967 *
968             I1 = N / 2
969             DO 120 J = 1, I1
970                LLWORK( J ) = .TRUE.
971   120       CONTINUE
972             DO 130 J = I1 + 1, N
973                LLWORK( J ) = .FALSE.
974   130       CONTINUE
975 *
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,
980      $            JTYPE, IOLDSD
981                INFO = ABS( IINFO )
982                GO TO 210
983             END IF
984 *
985             I1 = IN
986             DO 140 J = 1, I1
987                LLWORK( J ) = .FALSE.
988   140       CONTINUE
989             DO 150 J = I1 + 1, N
990                LLWORK( J ) = .TRUE.
991   150       CONTINUE
992 *
993             CALL DTGEVC( 'L', 'S', LLWORK, N, S1, LDA, P1, LDA,
994      $                   EVECTL( 1, I1+1 ), LDU, DUMMA, LDU, N, IN,
995      $                   WORK, IINFO )
996             IF( IINFO.NE.0 ) THEN
997                WRITE( NOUNIT, FMT = 9999 )'DTGEVC(L,S2)', IINFO, N,
998      $            JTYPE, IOLDSD
999                INFO = ABS( IINFO )
1000                GO TO 210
1001             END IF
1002 *
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
1009             END IF
1010 *
1011 *           10: Compute the left eigenvector Matrix with
1012 *               back transforming:
1013 *
1014             NTEST = 10
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,
1021      $            JTYPE, IOLDSD
1022                INFO = ABS( IINFO )
1023                GO TO 210
1024             END IF
1025 *
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
1032             END IF
1033 *
1034 *           11: Compute the right eigenvector Matrix without
1035 *               back transforming:
1036 *
1037             NTEST = 11
1038             RESULT( 11 ) = ULPINV
1039 *
1040 *           To test "SELECT" option, compute half of the eigenvectors
1041 *           in one call, and half in another
1042 *
1043             I1 = N / 2
1044             DO 160 J = 1, I1
1045                LLWORK( J ) = .TRUE.
1046   160       CONTINUE
1047             DO 170 J = I1 + 1, N
1048                LLWORK( J ) = .FALSE.
1049   170       CONTINUE
1050 *
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,
1055      $            JTYPE, IOLDSD
1056                INFO = ABS( IINFO )
1057                GO TO 210
1058             END IF
1059 *
1060             I1 = IN
1061             DO 180 J = 1, I1
1062                LLWORK( J ) = .FALSE.
1063   180       CONTINUE
1064             DO 190 J = I1 + 1, N
1065                LLWORK( J ) = .TRUE.
1066   190       CONTINUE
1067 *
1068             CALL DTGEVC( 'R', 'S', LLWORK, N, S1, LDA, P1, LDA, DUMMA,
1069      $                   LDU, EVECTR( 1, I1+1 ), LDU, N, IN, WORK,
1070      $                   IINFO )
1071             IF( IINFO.NE.0 ) THEN
1072                WRITE( NOUNIT, FMT = 9999 )'DTGEVC(R,S2)', IINFO, N,
1073      $            JTYPE, IOLDSD
1074                INFO = ABS( IINFO )
1075                GO TO 210
1076             END IF
1077 *
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
1084             END IF
1085 *
1086 *           12: Compute the right eigenvector Matrix with
1087 *               back transforming:
1088 *
1089             NTEST = 12
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,
1096      $            JTYPE, IOLDSD
1097                INFO = ABS( IINFO )
1098                GO TO 210
1099             END IF
1100 *
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
1107             END IF
1108 *
1109 *           Tests 13--15 are done only on request
1110 *
1111             IF( TSTDIF ) THEN
1112 *
1113 *              Do Tests 13--14
1114 *
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 ) )
1119 *
1120 *              Do Test 15
1121 *
1122                TEMP1 = ZERO
1123                TEMP2 = ZERO
1124                DO 200 J = 1, N
1125                   TEMP1 = MAX( TEMP1, ABS( ALPHR1( J )-ALPHR3( J ) )+
1126      $                    ABS( ALPHI1( J )-ALPHI3( J ) ) )
1127                   TEMP2 = MAX( TEMP2, ABS( BETA1( J )-BETA3( J ) ) )
1128   200          CONTINUE
1129 *
1130                TEMP1 = TEMP1 / MAX( SAFMIN, ULP*MAX( TEMP1, ANORM ) )
1131                TEMP2 = TEMP2 / MAX( SAFMIN, ULP*MAX( TEMP2, BNORM ) )
1132                RESULT( 15 ) = MAX( TEMP1, TEMP2 )
1133                NTEST = 15
1134             ELSE
1135                RESULT( 13 ) = ZERO
1136                RESULT( 14 ) = ZERO
1137                RESULT( 15 ) = ZERO
1138                NTEST = 12
1139             END IF
1140 *
1141 *           End of Loop -- Check for RESULT(j) > THRESH
1142 *
1143   210       CONTINUE
1144 *
1145             NTESTT = NTESTT + NTEST
1146 *
1147 *           Print out tests which fail.
1148 *
1149             DO 220 JR = 1, NTEST
1150                IF( RESULT( JR ).GE.THRESH ) THEN
1151 *
1152 *                 If this is the first test to fail,
1153 *                 print a header to the data file.
1154 *
1155                   IF( NERRS.EQ.0 ) THEN
1156                      WRITE( NOUNIT, FMT = 9997 )'DGG'
1157 *
1158 *                    Matrix types
1159 *
1160                      WRITE( NOUNIT, FMT = 9996 )
1161                      WRITE( NOUNIT, FMT = 9995 )
1162                      WRITE( NOUNIT, FMT = 9994 )'Orthogonal'
1163 *
1164 *                    Tests performed
1165 *
1166                      WRITE( NOUNIT, FMT = 9993 )'orthogonal', '''',
1167      $                  'transpose', ( '''', J = 1, 10 )
1168 *
1169                   END IF
1170                   NERRS = NERRS + 1
1171                   IF( RESULT( JR ).LT.10000.0D0 ) THEN
1172                      WRITE( NOUNIT, FMT = 9992 )N, JTYPE, IOLDSD, JR,
1173      $                  RESULT( JR )
1174                   ELSE
1175                      WRITE( NOUNIT, FMT = 9991 )N, JTYPE, IOLDSD, JR,
1176      $                  RESULT( JR )
1177                   END IF
1178                END IF
1179   220       CONTINUE
1180 *
1181   230    CONTINUE
1182   240 CONTINUE
1183 *
1184 *     Summary
1185 *
1186       CALL DLASUM( 'DGG', NOUNIT, NERRS, NTESTT )
1187       RETURN
1188 *
1189  9999 FORMAT( ' DCHKGG: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
1190      $      I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' )
1191 *
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,
1195      $      ')' )
1196 *
1197  9997 FORMAT( / 1X, A3, ' -- Real Generalized eigenvalue problem' )
1198 *
1199  9996 FORMAT( ' Matrix types (see DCHKGG for details): ' )
1200 *
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.' )
1217 *
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 )
1236 *
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 )
1241 *
1242 *     End of DCHKGG
1243 *
1244       END