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