1*> \brief \b DDRVGG
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 DDRVGG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12*                          THRSHN, NOUNIT, A, LDA, B, S, T, S2, T2, Q,
13*                          LDQ, Z, ALPHR1, ALPHI1, BETA1, ALPHR2, ALPHI2,
14*                          BETA2, VL, VR, WORK, LWORK, RESULT, INFO )
15*
16*       .. Scalar Arguments ..
17*       INTEGER            INFO, LDA, LDQ, LWORK, NOUNIT, NSIZES, NTYPES
18*       DOUBLE PRECISION   THRESH, THRSHN
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            ISEED( 4 ), NN( * )
23*       DOUBLE PRECISION   A( LDA, * ), ALPHI1( * ), ALPHI2( * ),
24*      $                   ALPHR1( * ), ALPHR2( * ), B( LDA, * ),
25*      $                   BETA1( * ), BETA2( * ), Q( LDQ, * ),
26*      $                   RESULT( * ), S( LDA, * ), S2( LDA, * ),
27*      $                   T( LDA, * ), T2( LDA, * ), VL( LDQ, * ),
28*      $                   VR( LDQ, * ), WORK( * ), Z( LDQ, * )
29*       ..
30*
31*
32*> \par Purpose:
33*  =============
34*>
35*> \verbatim
36*>
37*> DDRVGG  checks the nonsymmetric generalized eigenvalue driver
38*> routines.
39*>                               T          T        T
40*> DGEGS factors A and B as Q S Z  and Q T Z , where   means
41*> transpose, T is upper triangular, S is in generalized Schur form
42*> (block upper triangular, with 1x1 and 2x2 blocks on the diagonal,
43*> the 2x2 blocks corresponding to complex conjugate pairs of
44*> generalized eigenvalues), and Q and Z are orthogonal.  It also
45*> computes the generalized eigenvalues (alpha(1),beta(1)), ...,
46*> (alpha(n),beta(n)), where alpha(j)=S(j,j) and beta(j)=P(j,j) --
47*> thus, w(j) = alpha(j)/beta(j) is a root of the generalized
48*> eigenvalue problem
49*>
50*>     det( A - w(j) B ) = 0
51*>
52*> and m(j) = beta(j)/alpha(j) is a root of the essentially equivalent
53*> problem
54*>
55*>     det( m(j) A - B ) = 0
56*>
57*> DGEGV computes the generalized eigenvalues (alpha(1),beta(1)), ...,
58*> (alpha(n),beta(n)), the matrix L whose columns contain the
59*> generalized left eigenvectors l, and the matrix R whose columns
60*> contain the generalized right eigenvectors r for the pair (A,B).
61*>
62*> When DDRVGG is called, a number of matrix "sizes" ("n's") and a
63*> number of matrix "types" are specified.  For each size ("n")
64*> and each type of matrix, one matrix will be generated and used
65*> to test the nonsymmetric eigenroutines.  For each matrix, 7
66*> tests will be performed and compared with the threshhold THRESH:
67*>
68*> Results from DGEGS:
69*>
70*>                  T
71*> (1)   | A - Q S Z  | / ( |A| n ulp )
72*>
73*>                  T
74*> (2)   | B - Q T Z  | / ( |B| n ulp )
75*>
76*>               T
77*> (3)   | I - QQ  | / ( n ulp )
78*>
79*>               T
80*> (4)   | I - ZZ  | / ( n ulp )
81*>
82*> (5)   maximum over j of D(j)  where:
83*>
84*> if alpha(j) is real:
85*>                     |alpha(j) - S(j,j)|        |beta(j) - T(j,j)|
86*>           D(j) = ------------------------ + -----------------------
87*>                  max(|alpha(j)|,|S(j,j)|)   max(|beta(j)|,|T(j,j)|)
88*>
89*> if alpha(j) is complex:
90*>                                 | det( s S - w T ) |
91*>           D(j) = ---------------------------------------------------
92*>                  ulp max( s norm(S), |w| norm(T) )*norm( s S - w T )
93*>
94*>           and S and T are here the 2 x 2 diagonal blocks of S and T
95*>           corresponding to the j-th eigenvalue.
96*>
97*> Results from DGEGV:
98*>
99*> (6)   max over all left eigenvalue/-vector pairs (beta/alpha,l) of
100*>
101*>    | l**H * (beta A - alpha B) | / ( ulp max( |beta A|, |alpha B| ) )
102*>
103*>       where l**H is the conjugate tranpose of l.
104*>
105*> (7)   max over all right eigenvalue/-vector pairs (beta/alpha,r) of
106*>
107*>       | (beta A - alpha B) r | / ( ulp max( |beta A|, |alpha B| ) )
108*>
109*> Test Matrices
110*> ---- --------
111*>
112*> The sizes of the test matrices are specified by an array
113*> NN(1:NSIZES); the value of each element NN(j) specifies one size.
114*> The "types" are specified by a logical array DOTYPE( 1:NTYPES ); if
115*> DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
116*> Currently, the list of possible types is:
117*>
118*> (1)  ( 0, 0 )         (a pair of zero matrices)
119*>
120*> (2)  ( I, 0 )         (an identity and a zero matrix)
121*>
122*> (3)  ( 0, I )         (an identity and a zero matrix)
123*>
124*> (4)  ( I, I )         (a pair of identity matrices)
125*>
126*>         t   t
127*> (5)  ( J , J  )       (a pair of transposed Jordan blocks)
128*>
129*>                                     t                ( I   0  )
130*> (6)  ( X, Y )         where  X = ( J   0  )  and Y = (      t )
131*>                                  ( 0   I  )          ( 0   J  )
132*>                       and I is a k x k identity and J a (k+1)x(k+1)
133*>                       Jordan block; k=(N-1)/2
134*>
135*> (7)  ( D, I )         where D is diag( 0, 1,..., N-1 ) (a diagonal
136*>                       matrix with those diagonal entries.)
137*> (8)  ( I, D )
138*>
139*> (9)  ( big*D, small*I ) where "big" is near overflow and small=1/big
140*>
141*> (10) ( small*D, big*I )
142*>
143*> (11) ( big*I, small*D )
144*>
145*> (12) ( small*I, big*D )
146*>
147*> (13) ( big*D, big*I )
148*>
149*> (14) ( small*D, small*I )
150*>
151*> (15) ( D1, D2 )        where D1 is diag( 0, 0, 1, ..., N-3, 0 ) and
152*>                        D2 is diag( 0, N-3, N-4,..., 1, 0, 0 )
153*>           t   t
154*> (16) Q ( J , J ) Z     where Q and Z are random orthogonal matrices.
155*>
156*> (17) Q ( T1, T2 ) Z    where T1 and T2 are upper triangular matrices
157*>                        with random O(1) entries above the diagonal
158*>                        and diagonal entries diag(T1) =
159*>                        ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
160*>                        ( 0, N-3, N-4,..., 1, 0, 0 )
161*>
162*> (18) Q ( T1, T2 ) Z    diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
163*>                        diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
164*>                        s = machine precision.
165*>
166*> (19) Q ( T1, T2 ) Z    diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 )
167*>                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 )
168*>
169*>                                                        N-5
170*> (20) Q ( T1, T2 ) Z    diag(T1)=( 0, 0, 1, 1, a, ..., a   =s, 0 )
171*>                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
172*>
173*> (21) Q ( T1, T2 ) Z    diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 )
174*>                        diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
175*>                        where r1,..., r(N-4) are random.
176*>
177*> (22) Q ( big*T1, small*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
178*>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
179*>
180*> (23) Q ( small*T1, big*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
181*>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
182*>
183*> (24) Q ( small*T1, small*T2 ) Z  diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
184*>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
185*>
186*> (25) Q ( big*T1, big*T2 ) Z      diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
187*>                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
188*>
189*> (26) Q ( T1, T2 ) Z     where T1 and T2 are random upper-triangular
190*>                         matrices.
191*> \endverbatim
192*
193*  Arguments:
194*  ==========
195*
196*> \param[in] NSIZES
197*> \verbatim
198*>          NSIZES is INTEGER
199*>          The number of sizes of matrices to use.  If it is zero,
200*>          DDRVGG does nothing.  It must be at least zero.
201*> \endverbatim
202*>
203*> \param[in] NN
204*> \verbatim
205*>          NN is INTEGER array, dimension (NSIZES)
206*>          An array containing the sizes to be used for the matrices.
207*>          Zero values will be skipped.  The values must be at least
208*>          zero.
209*> \endverbatim
210*>
211*> \param[in] NTYPES
212*> \verbatim
213*>          NTYPES is INTEGER
214*>          The number of elements in DOTYPE.   If it is zero, DDRVGG
215*>          does nothing.  It must be at least zero.  If it is MAXTYP+1
216*>          and NSIZES is 1, then an additional type, MAXTYP+1 is
217*>          defined, which is to use whatever matrix is in A.  This
218*>          is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
219*>          DOTYPE(MAXTYP+1) is .TRUE. .
220*> \endverbatim
221*>
222*> \param[in] DOTYPE
223*> \verbatim
224*>          DOTYPE is LOGICAL array, dimension (NTYPES)
225*>          If DOTYPE(j) is .TRUE., then for each size in NN a
226*>          matrix of that size and of type j will be generated.
227*>          If NTYPES is smaller than the maximum number of types
228*>          defined (PARAMETER MAXTYP), then types NTYPES+1 through
229*>          MAXTYP will not be generated.  If NTYPES is larger
230*>          than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
231*>          will be ignored.
232*> \endverbatim
233*>
234*> \param[in,out] ISEED
235*> \verbatim
236*>          ISEED is INTEGER array, dimension (4)
237*>          On entry ISEED specifies the seed of the random number
238*>          generator. The array elements should be between 0 and 4095;
239*>          if not they will be reduced mod 4096.  Also, ISEED(4) must
240*>          be odd.  The random number generator uses a linear
241*>          congruential sequence limited to small integers, and so
242*>          should produce machine independent random numbers. The
243*>          values of ISEED are changed on exit, and can be used in the
244*>          next call to DDRVGG to continue the same random number
245*>          sequence.
246*> \endverbatim
247*>
248*> \param[in] THRESH
249*> \verbatim
250*>          THRESH is DOUBLE PRECISION
251*>          A test will count as "failed" if the "error", computed as
252*>          described above, exceeds THRESH.  Note that the error is
253*>          scaled to be O(1), so THRESH should be a reasonably small
254*>          multiple of 1, e.g., 10 or 100.  In particular, it should
255*>          not depend on the precision (single vs. double) or the size
256*>          of the matrix.  It must be at least zero.
257*> \endverbatim
258*>
259*> \param[in] THRSHN
260*> \verbatim
261*>          THRSHN is DOUBLE PRECISION
262*>          Threshhold for reporting eigenvector normalization error.
263*>          If the normalization of any eigenvector differs from 1 by
264*>          more than THRSHN*ulp, then a special error message will be
265*>          printed.  (This is handled separately from the other tests,
266*>          since only a compiler or programming error should cause an
267*>          error message, at least if THRSHN is at least 5--10.)
268*> \endverbatim
269*>
270*> \param[in] NOUNIT
271*> \verbatim
272*>          NOUNIT is INTEGER
273*>          The FORTRAN unit number for printing out error messages
274*>          (e.g., if a routine returns IINFO not equal to 0.)
275*> \endverbatim
276*>
277*> \param[in,out] A
278*> \verbatim
279*>          A is DOUBLE PRECISION array, dimension
280*>                            (LDA, max(NN))
281*>          Used to hold the original A matrix.  Used as input only
282*>          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
283*>          DOTYPE(MAXTYP+1)=.TRUE.
284*> \endverbatim
285*>
286*> \param[in] LDA
287*> \verbatim
288*>          LDA is INTEGER
289*>          The leading dimension of A, B, S, T, S2, and T2.
290*>          It must be at least 1 and at least max( NN ).
291*> \endverbatim
292*>
293*> \param[in,out] B
294*> \verbatim
295*>          B is DOUBLE PRECISION array, dimension
296*>                            (LDA, max(NN))
297*>          Used to hold the original B matrix.  Used as input only
298*>          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
299*>          DOTYPE(MAXTYP+1)=.TRUE.
300*> \endverbatim
301*>
302*> \param[out] S
303*> \verbatim
304*>          S is DOUBLE PRECISION array, dimension (LDA, max(NN))
305*>          The Schur form matrix computed from A by DGEGS.  On exit, S
306*>          contains the Schur form matrix corresponding to the matrix
307*>          in A.
308*> \endverbatim
309*>
310*> \param[out] T
311*> \verbatim
312*>          T is DOUBLE PRECISION array, dimension (LDA, max(NN))
313*>          The upper triangular matrix computed from B by DGEGS.
314*> \endverbatim
315*>
316*> \param[out] S2
317*> \verbatim
318*>          S2 is DOUBLE PRECISION array, dimension (LDA, max(NN))
319*>          The matrix computed from A by DGEGV.  This will be the
320*>          Schur form of some matrix related to A, but will not, in
321*>          general, be the same as S.
322*> \endverbatim
323*>
324*> \param[out] T2
325*> \verbatim
326*>          T2 is DOUBLE PRECISION array, dimension (LDA, max(NN))
327*>          The matrix computed from B by DGEGV.  This will be the
328*>          Schur form of some matrix related to B, but will not, in
329*>          general, be the same as T.
330*> \endverbatim
331*>
332*> \param[out] Q
333*> \verbatim
334*>          Q is DOUBLE PRECISION array, dimension (LDQ, max(NN))
335*>          The (left) orthogonal matrix computed by DGEGS.
336*> \endverbatim
337*>
338*> \param[in] LDQ
339*> \verbatim
340*>          LDQ is INTEGER
341*>          The leading dimension of Q, Z, VL, and VR.  It must
342*>          be at least 1 and at least max( NN ).
343*> \endverbatim
344*>
345*> \param[out] Z
346*> \verbatim
347*>          Z is DOUBLE PRECISION array of
348*>                             dimension( LDQ, max(NN) )
349*>          The (right) orthogonal matrix computed by DGEGS.
350*> \endverbatim
351*>
352*> \param[out] ALPHR1
353*> \verbatim
354*>          ALPHR1 is DOUBLE PRECISION array, dimension (max(NN))
355*> \endverbatim
356*>
357*> \param[out] ALPHI1
358*> \verbatim
359*>          ALPHI1 is DOUBLE PRECISION array, dimension (max(NN))
360*> \endverbatim
361*>
362*> \param[out] BETA1
363*> \verbatim
364*>          BETA1 is DOUBLE PRECISION array, dimension (max(NN))
365*>
366*>          The generalized eigenvalues of (A,B) computed by DGEGS.
367*>          ( ALPHR1(k)+ALPHI1(k)*i ) / BETA1(k) is the k-th
368*>          generalized eigenvalue of the matrices in A and B.
369*> \endverbatim
370*>
371*> \param[out] ALPHR2
372*> \verbatim
373*>          ALPHR2 is DOUBLE PRECISION array, dimension (max(NN))
374*> \endverbatim
375*>
376*> \param[out] ALPHI2
377*> \verbatim
378*>          ALPHI2 is DOUBLE PRECISION array, dimension (max(NN))
379*> \endverbatim
380*>
381*> \param[out] BETA2
382*> \verbatim
383*>          BETA2 is DOUBLE PRECISION array, dimension (max(NN))
384*>
385*>          The generalized eigenvalues of (A,B) computed by DGEGV.
386*>          ( ALPHR2(k)+ALPHI2(k)*i ) / BETA2(k) is the k-th
387*>          generalized eigenvalue of the matrices in A and B.
388*> \endverbatim
389*>
390*> \param[out] VL
391*> \verbatim
392*>          VL is DOUBLE PRECISION array, dimension (LDQ, max(NN))
393*>          The (block lower triangular) left eigenvector matrix for
394*>          the matrices in A and B.  (See DTGEVC for the format.)
395*> \endverbatim
396*>
397*> \param[out] VR
398*> \verbatim
399*>          VR is DOUBLE PRECISION array, dimension (LDQ, max(NN))
400*>          The (block upper triangular) right eigenvector matrix for
401*>          the matrices in A and B.  (See DTGEVC for the format.)
402*> \endverbatim
403*>
404*> \param[out] WORK
405*> \verbatim
406*>          WORK is DOUBLE PRECISION array, dimension (LWORK)
407*> \endverbatim
408*>
409*> \param[in] LWORK
410*> \verbatim
411*>          LWORK is INTEGER
412*>          The number of entries in WORK.  This must be at least
413*>          2*N + MAX( 6*N, N*(NB+1), (k+1)*(2*k+N+1) ), where
414*>          "k" is the sum of the blocksize and number-of-shifts for
415*>          DHGEQZ, and NB is the greatest of the blocksizes for
416*>          DGEQRF, DORMQR, and DORGQR.  (The blocksizes and the
417*>          number-of-shifts are retrieved through calls to ILAENV.)
418*> \endverbatim
419*>
420*> \param[out] RESULT
421*> \verbatim
422*>          RESULT is DOUBLE PRECISION array, dimension (15)
423*>          The values computed by the tests described above.
424*>          The values are currently limited to 1/ulp, to avoid
425*>          overflow.
426*> \endverbatim
427*>
428*> \param[out] INFO
429*> \verbatim
430*>          INFO is INTEGER
431*>          = 0:  successful exit
432*>          < 0:  if INFO = -i, the i-th argument had an illegal value.
433*>          > 0:  A routine returned an error code.  INFO is the
434*>                absolute value of the INFO value returned.
435*> \endverbatim
436*
437*  Authors:
438*  ========
439*
440*> \author Univ. of Tennessee
441*> \author Univ. of California Berkeley
442*> \author Univ. of Colorado Denver
443*> \author NAG Ltd.
444*
445*> \date November 2011
446*
447*> \ingroup double_eig
448*
449*  =====================================================================
450      SUBROUTINE DDRVGG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
451     $                   THRSHN, NOUNIT, A, LDA, B, S, T, S2, T2, Q,
452     $                   LDQ, Z, ALPHR1, ALPHI1, BETA1, ALPHR2, ALPHI2,
453     $                   BETA2, VL, VR, WORK, LWORK, RESULT, INFO )
454*
455*  -- LAPACK test routine (version 3.4.0) --
456*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
457*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
458*     November 2011
459*
460*     .. Scalar Arguments ..
461      INTEGER            INFO, LDA, LDQ, LWORK, NOUNIT, NSIZES, NTYPES
462      DOUBLE PRECISION   THRESH, THRSHN
463*     ..
464*     .. Array Arguments ..
465      LOGICAL            DOTYPE( * )
466      INTEGER            ISEED( 4 ), NN( * )
467      DOUBLE PRECISION   A( LDA, * ), ALPHI1( * ), ALPHI2( * ),
468     $                   ALPHR1( * ), ALPHR2( * ), B( LDA, * ),
469     $                   BETA1( * ), BETA2( * ), Q( LDQ, * ),
470     $                   RESULT( * ), S( LDA, * ), S2( LDA, * ),
471     $                   T( LDA, * ), T2( LDA, * ), VL( LDQ, * ),
472     $                   VR( LDQ, * ), WORK( * ), Z( LDQ, * )
473*     ..
474*
475*  =====================================================================
476*
477*     .. Parameters ..
478      DOUBLE PRECISION   ZERO, ONE
479      PARAMETER          ( ZERO = 0.0D0, ONE = 1.0D0 )
480      INTEGER            MAXTYP
481      PARAMETER          ( MAXTYP = 26 )
482*     ..
483*     .. Local Scalars ..
484      LOGICAL            BADNN, ILABAD
485      INTEGER            I1, IADD, IINFO, IN, J, JC, JR, JSIZE, JTYPE,
486     $                   LWKOPT, MTYPES, N, N1, NB, NBZ, NERRS, NMATS,
487     $                   NMAX, NS, NTEST, NTESTT
488      DOUBLE PRECISION   SAFMAX, SAFMIN, TEMP1, TEMP2, ULP, ULPINV
489*     ..
490*     .. Local Arrays ..
491      INTEGER            IASIGN( MAXTYP ), IBSIGN( MAXTYP ),
492     $                   IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ),
493     $                   KATYPE( MAXTYP ), KAZERO( MAXTYP ),
494     $                   KBMAGN( MAXTYP ), KBTYPE( MAXTYP ),
495     $                   KBZERO( MAXTYP ), KCLASS( MAXTYP ),
496     $                   KTRIAN( MAXTYP ), KZ1( 6 ), KZ2( 6 )
497      DOUBLE PRECISION   DUMMA( 4 ), RMAGN( 0: 3 )
498*     ..
499*     .. External Functions ..
500      INTEGER            ILAENV
501      DOUBLE PRECISION   DLAMCH, DLARND
502      EXTERNAL           ILAENV, DLAMCH, DLARND
503*     ..
504*     .. External Subroutines ..
505      EXTERNAL           ALASVM, DGEGS, DGEGV, DGET51, DGET52, DGET53,
506     $                   DLABAD, DLACPY, DLARFG, DLASET, DLATM4, DORM2R,
507     $                   XERBLA
508*     ..
509*     .. Intrinsic Functions ..
510      INTRINSIC          ABS, DBLE, MAX, MIN, SIGN
511*     ..
512*     .. Data statements ..
513      DATA               KCLASS / 15*1, 10*2, 1*3 /
514      DATA               KZ1 / 0, 1, 2, 1, 3, 3 /
515      DATA               KZ2 / 0, 0, 1, 2, 1, 1 /
516      DATA               KADD / 0, 0, 0, 0, 3, 2 /
517      DATA               KATYPE / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
518     $                   4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 /
519      DATA               KBTYPE / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
520     $                   1, 1, -4, 2, -4, 8*8, 0 /
521      DATA               KAZERO / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
522     $                   4*5, 4*3, 1 /
523      DATA               KBZERO / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
524     $                   4*6, 4*4, 1 /
525      DATA               KAMAGN / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
526     $                   2, 1 /
527      DATA               KBMAGN / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
528     $                   2, 1 /
529      DATA               KTRIAN / 16*0, 10*1 /
530      DATA               IASIGN / 6*0, 2, 0, 2*2, 2*0, 3*2, 0, 2, 3*0,
531     $                   5*2, 0 /
532      DATA               IBSIGN / 7*0, 2, 2*0, 2*2, 2*0, 2, 0, 2, 9*0 /
533*     ..
534*     .. Executable Statements ..
535*
536*     Check for errors
537*
538      INFO = 0
539*
540      BADNN = .FALSE.
541      NMAX = 1
542      DO 10 J = 1, NSIZES
543         NMAX = MAX( NMAX, NN( J ) )
544         IF( NN( J ).LT.0 )
545     $      BADNN = .TRUE.
546   10 CONTINUE
547*
548*     Maximum blocksize and shift -- we assume that blocksize and number
549*     of shifts are monotone increasing functions of N.
550*
551      NB = MAX( 1, ILAENV( 1, 'DGEQRF', ' ', NMAX, NMAX, -1, -1 ),
552     $     ILAENV( 1, 'DORMQR', 'LT', NMAX, NMAX, NMAX, -1 ),
553     $     ILAENV( 1, 'DORGQR', ' ', NMAX, NMAX, NMAX, -1 ) )
554      NBZ = ILAENV( 1, 'DHGEQZ', 'SII', NMAX, 1, NMAX, 0 )
555      NS = ILAENV( 4, 'DHGEQZ', 'SII', NMAX, 1, NMAX, 0 )
556      I1 = NBZ + NS
557      LWKOPT = 2*NMAX + MAX( 6*NMAX, NMAX*( NB+1 ),
558     $         ( 2*I1+NMAX+1 )*( I1+1 ) )
559*
560*     Check for errors
561*
562      IF( NSIZES.LT.0 ) THEN
563         INFO = -1
564      ELSE IF( BADNN ) THEN
565         INFO = -2
566      ELSE IF( NTYPES.LT.0 ) THEN
567         INFO = -3
568      ELSE IF( THRESH.LT.ZERO ) THEN
569         INFO = -6
570      ELSE IF( LDA.LE.1 .OR. LDA.LT.NMAX ) THEN
571         INFO = -10
572      ELSE IF( LDQ.LE.1 .OR. LDQ.LT.NMAX ) THEN
573         INFO = -19
574      ELSE IF( LWKOPT.GT.LWORK ) THEN
575         INFO = -30
576      END IF
577*
578      IF( INFO.NE.0 ) THEN
579         CALL XERBLA( 'DDRVGG', -INFO )
580         RETURN
581      END IF
582*
583*     Quick return if possible
584*
585      IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 )
586     $   RETURN
587*
588      SAFMIN = DLAMCH( 'Safe minimum' )
589      ULP = DLAMCH( 'Epsilon' )*DLAMCH( 'Base' )
590      SAFMIN = SAFMIN / ULP
591      SAFMAX = ONE / SAFMIN
592      CALL DLABAD( SAFMIN, SAFMAX )
593      ULPINV = ONE / ULP
594*
595*     The values RMAGN(2:3) depend on N, see below.
596*
597      RMAGN( 0 ) = ZERO
598      RMAGN( 1 ) = ONE
599*
600*     Loop over sizes, types
601*
602      NTESTT = 0
603      NERRS = 0
604      NMATS = 0
605*
606      DO 170 JSIZE = 1, NSIZES
607         N = NN( JSIZE )
608         N1 = MAX( 1, N )
609         RMAGN( 2 ) = SAFMAX*ULP / DBLE( N1 )
610         RMAGN( 3 ) = SAFMIN*ULPINV*N1
611*
612         IF( NSIZES.NE.1 ) THEN
613            MTYPES = MIN( MAXTYP, NTYPES )
614         ELSE
615            MTYPES = MIN( MAXTYP+1, NTYPES )
616         END IF
617*
618         DO 160 JTYPE = 1, MTYPES
619            IF( .NOT.DOTYPE( JTYPE ) )
620     $         GO TO 160
621            NMATS = NMATS + 1
622            NTEST = 0
623*
624*           Save ISEED in case of an error.
625*
626            DO 20 J = 1, 4
627               IOLDSD( J ) = ISEED( J )
628   20       CONTINUE
629*
630*           Initialize RESULT
631*
632            DO 30 J = 1, 15
633               RESULT( J ) = ZERO
634   30       CONTINUE
635*
636*           Compute A and B
637*
638*           Description of control parameters:
639*
640*           KZLASS: =1 means w/o rotation, =2 means w/ rotation,
641*                   =3 means random.
642*           KATYPE: the "type" to be passed to DLATM4 for computing A.
643*           KAZERO: the pattern of zeros on the diagonal for A:
644*                   =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
645*                   =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
646*                   =6: ( 0, 1, 0, xxx, 0 ).  (xxx means a string of
647*                   non-zero entries.)
648*           KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
649*                   =2: large, =3: small.
650*           IASIGN: 1 if the diagonal elements of A are to be
651*                   multiplied by a random magnitude 1 number, =2 if
652*                   randomly chosen diagonal blocks are to be rotated
653*                   to form 2x2 blocks.
654*           KBTYPE, KBZERO, KBMAGN, IBSIGN: the same, but for B.
655*           KTRIAN: =0: don't fill in the upper triangle, =1: do.
656*           KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
657*           RMAGN: used to implement KAMAGN and KBMAGN.
658*
659            IF( MTYPES.GT.MAXTYP )
660     $         GO TO 110
661            IINFO = 0
662            IF( KCLASS( JTYPE ).LT.3 ) THEN
663*
664*              Generate A (w/o rotation)
665*
666               IF( ABS( KATYPE( JTYPE ) ).EQ.3 ) THEN
667                  IN = 2*( ( N-1 ) / 2 ) + 1
668                  IF( IN.NE.N )
669     $               CALL DLASET( 'Full', N, N, ZERO, ZERO, A, LDA )
670               ELSE
671                  IN = N
672               END IF
673               CALL DLATM4( KATYPE( JTYPE ), IN, KZ1( KAZERO( JTYPE ) ),
674     $                      KZ2( KAZERO( JTYPE ) ), IASIGN( JTYPE ),
675     $                      RMAGN( KAMAGN( JTYPE ) ), ULP,
676     $                      RMAGN( KTRIAN( JTYPE )*KAMAGN( JTYPE ) ), 2,
677     $                      ISEED, A, LDA )
678               IADD = KADD( KAZERO( JTYPE ) )
679               IF( IADD.GT.0 .AND. IADD.LE.N )
680     $            A( IADD, IADD ) = ONE
681*
682*              Generate B (w/o rotation)
683*
684               IF( ABS( KBTYPE( JTYPE ) ).EQ.3 ) THEN
685                  IN = 2*( ( N-1 ) / 2 ) + 1
686                  IF( IN.NE.N )
687     $               CALL DLASET( 'Full', N, N, ZERO, ZERO, B, LDA )
688               ELSE
689                  IN = N
690               END IF
691               CALL DLATM4( KBTYPE( JTYPE ), IN, KZ1( KBZERO( JTYPE ) ),
692     $                      KZ2( KBZERO( JTYPE ) ), IBSIGN( JTYPE ),
693     $                      RMAGN( KBMAGN( JTYPE ) ), ONE,
694     $                      RMAGN( KTRIAN( JTYPE )*KBMAGN( JTYPE ) ), 2,
695     $                      ISEED, B, LDA )
696               IADD = KADD( KBZERO( JTYPE ) )
697               IF( IADD.NE.0 .AND. IADD.LE.N )
698     $            B( IADD, IADD ) = ONE
699*
700               IF( KCLASS( JTYPE ).EQ.2 .AND. N.GT.0 ) THEN
701*
702*                 Include rotations
703*
704*                 Generate Q, Z as Householder transformations times
705*                 a diagonal matrix.
706*
707                  DO 50 JC = 1, N - 1
708                     DO 40 JR = JC, N
709                        Q( JR, JC ) = DLARND( 3, ISEED )
710                        Z( JR, JC ) = DLARND( 3, ISEED )
711   40                CONTINUE
712                     CALL DLARFG( N+1-JC, Q( JC, JC ), Q( JC+1, JC ), 1,
713     $                            WORK( JC ) )
714                     WORK( 2*N+JC ) = SIGN( ONE, Q( JC, JC ) )
715                     Q( JC, JC ) = ONE
716                     CALL DLARFG( N+1-JC, Z( JC, JC ), Z( JC+1, JC ), 1,
717     $                            WORK( N+JC ) )
718                     WORK( 3*N+JC ) = SIGN( ONE, Z( JC, JC ) )
719                     Z( JC, JC ) = ONE
720   50             CONTINUE
721                  Q( N, N ) = ONE
722                  WORK( N ) = ZERO
723                  WORK( 3*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
724                  Z( N, N ) = ONE
725                  WORK( 2*N ) = ZERO
726                  WORK( 4*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
727*
728*                 Apply the diagonal matrices
729*
730                  DO 70 JC = 1, N
731                     DO 60 JR = 1, N
732                        A( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
733     $                                A( JR, JC )
734                        B( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
735     $                                B( JR, JC )
736   60                CONTINUE
737   70             CONTINUE
738                  CALL DORM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, A,
739     $                         LDA, WORK( 2*N+1 ), IINFO )
740                  IF( IINFO.NE.0 )
741     $               GO TO 100
742                  CALL DORM2R( 'R', 'T', N, N, N-1, Z, LDQ, WORK( N+1 ),
743     $                         A, LDA, WORK( 2*N+1 ), IINFO )
744                  IF( IINFO.NE.0 )
745     $               GO TO 100
746                  CALL DORM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, B,
747     $                         LDA, WORK( 2*N+1 ), IINFO )
748                  IF( IINFO.NE.0 )
749     $               GO TO 100
750                  CALL DORM2R( 'R', 'T', N, N, N-1, Z, LDQ, WORK( N+1 ),
751     $                         B, LDA, WORK( 2*N+1 ), IINFO )
752                  IF( IINFO.NE.0 )
753     $               GO TO 100
754               END IF
755            ELSE
756*
757*              Random matrices
758*
759               DO 90 JC = 1, N
760                  DO 80 JR = 1, N
761                     A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
762     $                             DLARND( 2, ISEED )
763                     B( JR, JC ) = RMAGN( KBMAGN( JTYPE ) )*
764     $                             DLARND( 2, ISEED )
765   80             CONTINUE
766   90          CONTINUE
767            END IF
768*
769  100       CONTINUE
770*
771            IF( IINFO.NE.0 ) THEN
772               WRITE( NOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE,
773     $            IOLDSD
774               INFO = ABS( IINFO )
775               RETURN
776            END IF
777*
778  110       CONTINUE
779*
780*           Call DGEGS to compute H, T, Q, Z, alpha, and beta.
781*
782            CALL DLACPY( ' ', N, N, A, LDA, S, LDA )
783            CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
784            NTEST = 1
785            RESULT( 1 ) = ULPINV
786*
787            CALL DGEGS( 'V', 'V', N, S, LDA, T, LDA, ALPHR1, ALPHI1,
788     $                  BETA1, Q, LDQ, Z, LDQ, WORK, LWORK, IINFO )
789            IF( IINFO.NE.0 ) THEN
790               WRITE( NOUNIT, FMT = 9999 )'DGEGS', IINFO, N, JTYPE,
791     $            IOLDSD
792               INFO = ABS( IINFO )
793               GO TO 140
794            END IF
795*
796            NTEST = 4
797*
798*           Do tests 1--4
799*
800            CALL DGET51( 1, N, A, LDA, S, LDA, Q, LDQ, Z, LDQ, WORK,
801     $                   RESULT( 1 ) )
802            CALL DGET51( 1, N, B, LDA, T, LDA, Q, LDQ, Z, LDQ, WORK,
803     $                   RESULT( 2 ) )
804            CALL DGET51( 3, N, B, LDA, T, LDA, Q, LDQ, Q, LDQ, WORK,
805     $                   RESULT( 3 ) )
806            CALL DGET51( 3, N, B, LDA, T, LDA, Z, LDQ, Z, LDQ, WORK,
807     $                   RESULT( 4 ) )
808*
809*           Do test 5: compare eigenvalues with diagonals.
810*           Also check Schur form of A.
811*
812            TEMP1 = ZERO
813*
814            DO 120 J = 1, N
815               ILABAD = .FALSE.
816               IF( ALPHI1( J ).EQ.ZERO ) THEN
817                  TEMP2 = ( ABS( ALPHR1( J )-S( J, J ) ) /
818     $                    MAX( SAFMIN, ABS( ALPHR1( J ) ), ABS( S( J,
819     $                    J ) ) )+ABS( BETA1( J )-T( J, J ) ) /
820     $                    MAX( SAFMIN, ABS( BETA1( J ) ), ABS( T( J,
821     $                    J ) ) ) ) / ULP
822                  IF( J.LT.N ) THEN
823                     IF( S( J+1, J ).NE.ZERO )
824     $                  ILABAD = .TRUE.
825                  END IF
826                  IF( J.GT.1 ) THEN
827                     IF( S( J, J-1 ).NE.ZERO )
828     $                  ILABAD = .TRUE.
829                  END IF
830               ELSE
831                  IF( ALPHI1( J ).GT.ZERO ) THEN
832                     I1 = J
833                  ELSE
834                     I1 = J - 1
835                  END IF
836                  IF( I1.LE.0 .OR. I1.GE.N ) THEN
837                     ILABAD = .TRUE.
838                  ELSE IF( I1.LT.N-1 ) THEN
839                     IF( S( I1+2, I1+1 ).NE.ZERO )
840     $                  ILABAD = .TRUE.
841                  ELSE IF( I1.GT.1 ) THEN
842                     IF( S( I1, I1-1 ).NE.ZERO )
843     $                  ILABAD = .TRUE.
844                  END IF
845                  IF( .NOT.ILABAD ) THEN
846                     CALL DGET53( S( I1, I1 ), LDA, T( I1, I1 ), LDA,
847     $                            BETA1( J ), ALPHR1( J ), ALPHI1( J ),
848     $                            TEMP2, IINFO )
849                     IF( IINFO.GE.3 ) THEN
850                        WRITE( NOUNIT, FMT = 9997 )IINFO, J, N, JTYPE,
851     $                     IOLDSD
852                        INFO = ABS( IINFO )
853                     END IF
854                  ELSE
855                     TEMP2 = ULPINV
856                  END IF
857               END IF
858               TEMP1 = MAX( TEMP1, TEMP2 )
859               IF( ILABAD ) THEN
860                  WRITE( NOUNIT, FMT = 9996 )J, N, JTYPE, IOLDSD
861               END IF
862  120       CONTINUE
863            RESULT( 5 ) = TEMP1
864*
865*           Call DGEGV to compute S2, T2, VL, and VR, do tests.
866*
867*           Eigenvalues and Eigenvectors
868*
869            CALL DLACPY( ' ', N, N, A, LDA, S2, LDA )
870            CALL DLACPY( ' ', N, N, B, LDA, T2, LDA )
871            NTEST = 6
872            RESULT( 6 ) = ULPINV
873*
874            CALL DGEGV( 'V', 'V', N, S2, LDA, T2, LDA, ALPHR2, ALPHI2,
875     $                  BETA2, VL, LDQ, VR, LDQ, WORK, LWORK, IINFO )
876            IF( IINFO.NE.0 ) THEN
877               WRITE( NOUNIT, FMT = 9999 )'DGEGV', IINFO, N, JTYPE,
878     $            IOLDSD
879               INFO = ABS( IINFO )
880               GO TO 140
881            END IF
882*
883            NTEST = 7
884*
885*           Do Tests 6 and 7
886*
887            CALL DGET52( .TRUE., N, A, LDA, B, LDA, VL, LDQ, ALPHR2,
888     $                   ALPHI2, BETA2, WORK, DUMMA( 1 ) )
889            RESULT( 6 ) = DUMMA( 1 )
890            IF( DUMMA( 2 ).GT.THRSHN ) THEN
891               WRITE( NOUNIT, FMT = 9998 )'Left', 'DGEGV', DUMMA( 2 ),
892     $            N, JTYPE, IOLDSD
893            END IF
894*
895            CALL DGET52( .FALSE., N, A, LDA, B, LDA, VR, LDQ, ALPHR2,
896     $                   ALPHI2, BETA2, WORK, DUMMA( 1 ) )
897            RESULT( 7 ) = DUMMA( 1 )
898            IF( DUMMA( 2 ).GT.THRESH ) THEN
899               WRITE( NOUNIT, FMT = 9998 )'Right', 'DGEGV', DUMMA( 2 ),
900     $            N, JTYPE, IOLDSD
901            END IF
902*
903*           Check form of Complex eigenvalues.
904*
905            DO 130 J = 1, N
906               ILABAD = .FALSE.
907               IF( ALPHI2( J ).GT.ZERO ) THEN
908                  IF( J.EQ.N ) THEN
909                     ILABAD = .TRUE.
910                  ELSE IF( ALPHI2( J+1 ).GE.ZERO ) THEN
911                     ILABAD = .TRUE.
912                  END IF
913               ELSE IF( ALPHI2( J ).LT.ZERO ) THEN
914                  IF( J.EQ.1 ) THEN
915                     ILABAD = .TRUE.
916                  ELSE IF( ALPHI2( J-1 ).LE.ZERO ) THEN
917                     ILABAD = .TRUE.
918                  END IF
919               END IF
920               IF( ILABAD ) THEN
921                  WRITE( NOUNIT, FMT = 9996 )J, N, JTYPE, IOLDSD
922               END IF
923  130       CONTINUE
924*
925*           End of Loop -- Check for RESULT(j) > THRESH
926*
927  140       CONTINUE
928*
929            NTESTT = NTESTT + NTEST
930*
931*           Print out tests which fail.
932*
933            DO 150 JR = 1, NTEST
934               IF( RESULT( JR ).GE.THRESH ) THEN
935*
936*                 If this is the first test to fail,
937*                 print a header to the data file.
938*
939                  IF( NERRS.EQ.0 ) THEN
940                     WRITE( NOUNIT, FMT = 9995 )'DGG'
941*
942*                    Matrix types
943*
944                     WRITE( NOUNIT, FMT = 9994 )
945                     WRITE( NOUNIT, FMT = 9993 )
946                     WRITE( NOUNIT, FMT = 9992 )'Orthogonal'
947*
948*                    Tests performed
949*
950                     WRITE( NOUNIT, FMT = 9991 )'orthogonal', '''',
951     $                  'transpose', ( '''', J = 1, 5 )
952*
953                  END IF
954                  NERRS = NERRS + 1
955                  IF( RESULT( JR ).LT.10000.0D0 ) THEN
956                     WRITE( NOUNIT, FMT = 9990 )N, JTYPE, IOLDSD, JR,
957     $                  RESULT( JR )
958                  ELSE
959                     WRITE( NOUNIT, FMT = 9989 )N, JTYPE, IOLDSD, JR,
960     $                  RESULT( JR )
961                  END IF
962               END IF
963  150       CONTINUE
964*
965  160    CONTINUE
966  170 CONTINUE
967*
968*     Summary
969*
970      CALL ALASVM( 'DGG', NOUNIT, NERRS, NTESTT, 0 )
971      RETURN
972*
973 9999 FORMAT( ' DDRVGG: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
974     $      I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' )
975*
976 9998 FORMAT( ' DDRVGG: ', A, ' Eigenvectors from ', A, ' incorrectly ',
977     $      'normalized.', / ' Bits of error=', 0P, G10.3, ',', 9X,
978     $      'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5,
979     $      ')' )
980*
981 9997 FORMAT( ' DDRVGG: DGET53 returned INFO=', I1, ' for eigenvalue ',
982     $      I6, '.', / 9X, 'N=', I6, ', JTYPE=', I6, ', ISEED=(',
983     $      3( I5, ',' ), I5, ')' )
984*
985 9996 FORMAT( ' DDRVGG: S not in Schur form at eigenvalue ', I6, '.',
986     $      / 9X, 'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ),
987     $      I5, ')' )
988*
989 9995 FORMAT( / 1X, A3, ' -- Real Generalized eigenvalue problem driver'
990     $       )
991*
992 9994 FORMAT( ' Matrix types (see DDRVGG for details): ' )
993*
994 9993 FORMAT( ' Special Matrices:', 23X,
995     $      '(J''=transposed Jordan block)',
996     $      / '   1=(0,0)  2=(I,0)  3=(0,I)  4=(I,I)  5=(J'',J'')  ',
997     $      '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices:  ( ',
998     $      'D=diag(0,1,2,...) )', / '   7=(D,I)   9=(large*D, small*I',
999     $      ')  11=(large*I, small*D)  13=(large*D, large*I)', /
1000     $      '   8=(I,D)  10=(small*D, large*I)  12=(small*I, large*D) ',
1001     $      ' 14=(small*D, small*I)', / '  15=(D, reversed D)' )
1002 9992 FORMAT( ' Matrices Rotated by Random ', A, ' Matrices U, V:',
1003     $      / '  16=Transposed Jordan Blocks             19=geometric ',
1004     $      'alpha, beta=0,1', / '  17=arithm. alpha&beta             ',
1005     $      '      20=arithmetic alpha, beta=0,1', / '  18=clustered ',
1006     $      'alpha, beta=0,1            21=random alpha, beta=0,1',
1007     $      / ' Large & Small Matrices:', / '  22=(large, small)   ',
1008     $      '23=(small,large)    24=(small,small)    25=(large,large)',
1009     $      / '  26=random O(1) matrices.' )
1010*
1011 9991 FORMAT( / ' Tests performed:  (S is Schur, T is triangular, ',
1012     $      'Q and Z are ', A, ',', / 20X,
1013     $      'l and r are the appropriate left and right', / 19X,
1014     $      'eigenvectors, resp., a is alpha, b is beta, and', / 19X, A,
1015     $      ' means ', A, '.)', / ' 1 = | A - Q S Z', A,
1016     $      ' | / ( |A| n ulp )      2 = | B - Q T Z', A,
1017     $      ' | / ( |B| n ulp )', / ' 3 = | I - QQ', A,
1018     $      ' | / ( n ulp )             4 = | I - ZZ', A,
1019     $      ' | / ( n ulp )', /
1020     $      ' 5 = difference between (alpha,beta) and diagonals of',
1021     $      ' (S,T)', / ' 6 = max | ( b A - a B )', A,
1022     $      ' l | / const.   7 = max | ( b A - a B ) r | / const.',
1023     $      / 1X )
1024 9990 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
1025     $      4( I4, ',' ), ' result ', I3, ' is', 0P, F8.2 )
1026 9989 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
1027     $      4( I4, ',' ), ' result ', I3, ' is', 1P, D10.3 )
1028*
1029*     End of DDRVGG
1030*
1031      END
1032