LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchkhs.f
Go to the documentation of this file.
1*> \brief \b ZCHKHS
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 ZCHKHS( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12* NOUNIT, A, LDA, H, T1, T2, U, LDU, Z, UZ, W1,
13* W3, EVECTL, EVECTR, EVECTY, EVECTX, UU, TAU,
14* WORK, NWORK, RWORK, IWORK, SELECT, RESULT,
15* INFO )
16*
17* .. Scalar Arguments ..
18* INTEGER INFO, LDA, LDU, NOUNIT, NSIZES, NTYPES, NWORK
19* DOUBLE PRECISION THRESH
20* ..
21* .. Array Arguments ..
22* LOGICAL DOTYPE( * ), SELECT( * )
23* INTEGER ISEED( 4 ), IWORK( * ), NN( * )
24* DOUBLE PRECISION RESULT( 16 ), RWORK( * )
25* COMPLEX*16 A( LDA, * ), EVECTL( LDU, * ),
26* $ EVECTR( LDU, * ), EVECTX( LDU, * ),
27* $ EVECTY( LDU, * ), H( LDA, * ), T1( LDA, * ),
28* $ T2( LDA, * ), TAU( * ), U( LDU, * ),
29* $ UU( LDU, * ), UZ( LDU, * ), W1( * ), W3( * ),
30* $ WORK( * ), Z( LDU, * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> ZCHKHS checks the nonsymmetric eigenvalue problem routines.
40*>
41*> ZGEHRD factors A as U H U' , where ' means conjugate
42*> transpose, H is hessenberg, and U is unitary.
43*>
44*> ZUNGHR generates the unitary matrix U.
45*>
46*> ZUNMHR multiplies a matrix by the unitary matrix U.
47*>
48*> ZHSEQR factors H as Z T Z' , where Z is unitary and T
49*> is upper triangular. It also computes the eigenvalues,
50*> w(1), ..., w(n); we define a diagonal matrix W whose
51*> (diagonal) entries are the eigenvalues.
52*>
53*> ZTREVC computes the left eigenvector matrix L and the
54*> right eigenvector matrix R for the matrix T. The
55*> columns of L are the complex conjugates of the left
56*> eigenvectors of T. The columns of R are the right
57*> eigenvectors of T. L is lower triangular, and R is
58*> upper triangular.
59*>
60*> ZHSEIN computes the left eigenvector matrix Y and the
61*> right eigenvector matrix X for the matrix H. The
62*> columns of Y are the complex conjugates of the left
63*> eigenvectors of H. The columns of X are the right
64*> eigenvectors of H. Y is lower triangular, and X is
65*> upper triangular.
66*>
67*> ZTREVC3 computes left and right eigenvector matrices
68*> from a Schur matrix T and backtransforms them with Z
69*> to eigenvector matrices L and R for A. L and R are
70*> GE matrices.
71*>
72*> When ZCHKHS is called, a number of matrix "sizes" ("n's") and a
73*> number of matrix "types" are specified. For each size ("n")
74*> and each type of matrix, one matrix will be generated and used
75*> to test the nonsymmetric eigenroutines. For each matrix, 16
76*> tests will be performed:
77*>
78*> (1) | A - U H U**H | / ( |A| n ulp )
79*>
80*> (2) | I - UU**H | / ( n ulp )
81*>
82*> (3) | H - Z T Z**H | / ( |H| n ulp )
83*>
84*> (4) | I - ZZ**H | / ( n ulp )
85*>
86*> (5) | A - UZ H (UZ)**H | / ( |A| n ulp )
87*>
88*> (6) | I - UZ (UZ)**H | / ( n ulp )
89*>
90*> (7) | T(Z computed) - T(Z not computed) | / ( |T| ulp )
91*>
92*> (8) | W(Z computed) - W(Z not computed) | / ( |W| ulp )
93*>
94*> (9) | TR - RW | / ( |T| |R| ulp )
95*>
96*> (10) | L**H T - W**H L | / ( |T| |L| ulp )
97*>
98*> (11) | HX - XW | / ( |H| |X| ulp )
99*>
100*> (12) | Y**H H - W**H Y | / ( |H| |Y| ulp )
101*>
102*> (13) | AX - XW | / ( |A| |X| ulp )
103*>
104*> (14) | Y**H A - W**H Y | / ( |A| |Y| ulp )
105*>
106*> (15) | AR - RW | / ( |A| |R| ulp )
107*>
108*> (16) | LA - WL | / ( |A| |L| ulp )
109*>
110*> The "sizes" are specified by an array NN(1:NSIZES); the value of
111*> each element NN(j) specifies one size.
112*> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
113*> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
114*> Currently, the list of possible types is:
115*>
116*> (1) The zero matrix.
117*> (2) The identity matrix.
118*> (3) A (transposed) Jordan block, with 1's on the diagonal.
119*>
120*> (4) A diagonal matrix with evenly spaced entries
121*> 1, ..., ULP and random complex angles.
122*> (ULP = (first number larger than 1) - 1 )
123*> (5) A diagonal matrix with geometrically spaced entries
124*> 1, ..., ULP and random complex angles.
125*> (6) A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
126*> and random complex angles.
127*>
128*> (7) Same as (4), but multiplied by SQRT( overflow threshold )
129*> (8) Same as (4), but multiplied by SQRT( underflow threshold )
130*>
131*> (9) A matrix of the form U' T U, where U is unitary and
132*> T has evenly spaced entries 1, ..., ULP with random complex
133*> angles on the diagonal and random O(1) entries in the upper
134*> triangle.
135*>
136*> (10) A matrix of the form U' T U, where U is unitary and
137*> T has geometrically spaced entries 1, ..., ULP with random
138*> complex angles on the diagonal and random O(1) entries in
139*> the upper triangle.
140*>
141*> (11) A matrix of the form U' T U, where U is unitary and
142*> T has "clustered" entries 1, ULP,..., ULP with random
143*> complex angles on the diagonal and random O(1) entries in
144*> the upper triangle.
145*>
146*> (12) A matrix of the form U' T U, where U is unitary and
147*> T has complex eigenvalues randomly chosen from
148*> ULP < |z| < 1 and random O(1) entries in the upper
149*> triangle.
150*>
151*> (13) A matrix of the form X' T X, where X has condition
152*> SQRT( ULP ) and T has evenly spaced entries 1, ..., ULP
153*> with random complex angles on the diagonal and random O(1)
154*> entries in the upper triangle.
155*>
156*> (14) A matrix of the form X' T X, where X has condition
157*> SQRT( ULP ) and T has geometrically spaced entries
158*> 1, ..., ULP with random complex angles on the diagonal
159*> and random O(1) entries in the upper triangle.
160*>
161*> (15) A matrix of the form X' T X, where X has condition
162*> SQRT( ULP ) and T has "clustered" entries 1, ULP,..., ULP
163*> with random complex angles on the diagonal and random O(1)
164*> entries in the upper triangle.
165*>
166*> (16) A matrix of the form X' T X, where X has condition
167*> SQRT( ULP ) and T has complex eigenvalues randomly chosen
168*> from ULP < |z| < 1 and random O(1) entries in the upper
169*> triangle.
170*>
171*> (17) Same as (16), but multiplied by SQRT( overflow threshold )
172*> (18) Same as (16), but multiplied by SQRT( underflow threshold )
173*>
174*> (19) Nonsymmetric matrix with random entries chosen from |z| < 1
175*> (20) Same as (19), but multiplied by SQRT( overflow threshold )
176*> (21) Same as (19), but multiplied by SQRT( underflow threshold )
177*> \endverbatim
178*
179* Arguments:
180* ==========
181*
182*> \verbatim
183*> NSIZES - INTEGER
184*> The number of sizes of matrices to use. If it is zero,
185*> ZCHKHS does nothing. It must be at least zero.
186*> Not modified.
187*>
188*> NN - INTEGER array, dimension (NSIZES)
189*> An array containing the sizes to be used for the matrices.
190*> Zero values will be skipped. The values must be at least
191*> zero.
192*> Not modified.
193*>
194*> NTYPES - INTEGER
195*> The number of elements in DOTYPE. If it is zero, ZCHKHS
196*> does nothing. It must be at least zero. If it is MAXTYP+1
197*> and NSIZES is 1, then an additional type, MAXTYP+1 is
198*> defined, which is to use whatever matrix is in A. This
199*> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
200*> DOTYPE(MAXTYP+1) is .TRUE. .
201*> Not modified.
202*>
203*> DOTYPE - LOGICAL array, dimension (NTYPES)
204*> If DOTYPE(j) is .TRUE., then for each size in NN a
205*> matrix of that size and of type j will be generated.
206*> If NTYPES is smaller than the maximum number of types
207*> defined (PARAMETER MAXTYP), then types NTYPES+1 through
208*> MAXTYP will not be generated. If NTYPES is larger
209*> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
210*> will be ignored.
211*> Not modified.
212*>
213*> ISEED - INTEGER array, dimension (4)
214*> On entry ISEED specifies the seed of the random number
215*> generator. The array elements should be between 0 and 4095;
216*> if not they will be reduced mod 4096. Also, ISEED(4) must
217*> be odd. The random number generator uses a linear
218*> congruential sequence limited to small integers, and so
219*> should produce machine independent random numbers. The
220*> values of ISEED are changed on exit, and can be used in the
221*> next call to ZCHKHS to continue the same random number
222*> sequence.
223*> Modified.
224*>
225*> THRESH - DOUBLE PRECISION
226*> A test will count as "failed" if the "error", computed as
227*> described above, exceeds THRESH. Note that the error
228*> is scaled to be O(1), so THRESH should be a reasonably
229*> small multiple of 1, e.g., 10 or 100. In particular,
230*> it should not depend on the precision (single vs. double)
231*> or the size of the matrix. It must be at least zero.
232*> Not modified.
233*>
234*> NOUNIT - INTEGER
235*> The FORTRAN unit number for printing out error messages
236*> (e.g., if a routine returns IINFO not equal to 0.)
237*> Not modified.
238*>
239*> A - COMPLEX*16 array, dimension (LDA,max(NN))
240*> Used to hold the matrix whose eigenvalues are to be
241*> computed. On exit, A contains the last matrix actually
242*> used.
243*> Modified.
244*>
245*> LDA - INTEGER
246*> The leading dimension of A, H, T1 and T2. It must be at
247*> least 1 and at least max( NN ).
248*> Not modified.
249*>
250*> H - COMPLEX*16 array, dimension (LDA,max(NN))
251*> The upper hessenberg matrix computed by ZGEHRD. On exit,
252*> H contains the Hessenberg form of the matrix in A.
253*> Modified.
254*>
255*> T1 - COMPLEX*16 array, dimension (LDA,max(NN))
256*> The Schur (="quasi-triangular") matrix computed by ZHSEQR
257*> if Z is computed. On exit, T1 contains the Schur form of
258*> the matrix in A.
259*> Modified.
260*>
261*> T2 - COMPLEX*16 array, dimension (LDA,max(NN))
262*> The Schur matrix computed by ZHSEQR when Z is not computed.
263*> This should be identical to T1.
264*> Modified.
265*>
266*> LDU - INTEGER
267*> The leading dimension of U, Z, UZ and UU. It must be at
268*> least 1 and at least max( NN ).
269*> Not modified.
270*>
271*> U - COMPLEX*16 array, dimension (LDU,max(NN))
272*> The unitary matrix computed by ZGEHRD.
273*> Modified.
274*>
275*> Z - COMPLEX*16 array, dimension (LDU,max(NN))
276*> The unitary matrix computed by ZHSEQR.
277*> Modified.
278*>
279*> UZ - COMPLEX*16 array, dimension (LDU,max(NN))
280*> The product of U times Z.
281*> Modified.
282*>
283*> W1 - COMPLEX*16 array, dimension (max(NN))
284*> The eigenvalues of A, as computed by a full Schur
285*> decomposition H = Z T Z'. On exit, W1 contains the
286*> eigenvalues of the matrix in A.
287*> Modified.
288*>
289*> W3 - COMPLEX*16 array, dimension (max(NN))
290*> The eigenvalues of A, as computed by a partial Schur
291*> decomposition (Z not computed, T only computed as much
292*> as is necessary for determining eigenvalues). On exit,
293*> W3 contains the eigenvalues of the matrix in A, possibly
294*> perturbed by ZHSEIN.
295*> Modified.
296*>
297*> EVECTL - COMPLEX*16 array, dimension (LDU,max(NN))
298*> The conjugate transpose of the (upper triangular) left
299*> eigenvector matrix for the matrix in T1.
300*> Modified.
301*>
302*> EVEZTR - COMPLEX*16 array, dimension (LDU,max(NN))
303*> The (upper triangular) right eigenvector matrix for the
304*> matrix in T1.
305*> Modified.
306*>
307*> EVECTY - COMPLEX*16 array, dimension (LDU,max(NN))
308*> The conjugate transpose of the left eigenvector matrix
309*> for the matrix in H.
310*> Modified.
311*>
312*> EVECTX - COMPLEX*16 array, dimension (LDU,max(NN))
313*> The right eigenvector matrix for the matrix in H.
314*> Modified.
315*>
316*> UU - COMPLEX*16 array, dimension (LDU,max(NN))
317*> Details of the unitary matrix computed by ZGEHRD.
318*> Modified.
319*>
320*> TAU - COMPLEX*16 array, dimension (max(NN))
321*> Further details of the unitary matrix computed by ZGEHRD.
322*> Modified.
323*>
324*> WORK - COMPLEX*16 array, dimension (NWORK)
325*> Workspace.
326*> Modified.
327*>
328*> NWORK - INTEGER
329*> The number of entries in WORK. NWORK >= 4*NN(j)*NN(j) + 2.
330*>
331*> RWORK - DOUBLE PRECISION array, dimension (max(NN))
332*> Workspace. Could be equivalenced to IWORK, but not SELECT.
333*> Modified.
334*>
335*> IWORK - INTEGER array, dimension (max(NN))
336*> Workspace.
337*> Modified.
338*>
339*> SELECT - LOGICAL array, dimension (max(NN))
340*> Workspace. Could be equivalenced to IWORK, but not RWORK.
341*> Modified.
342*>
343*> RESULT - DOUBLE PRECISION array, dimension (16)
344*> The values computed by the fourteen tests described above.
345*> The values are currently limited to 1/ulp, to avoid
346*> overflow.
347*> Modified.
348*>
349*> INFO - INTEGER
350*> If 0, then everything ran OK.
351*> -1: NSIZES < 0
352*> -2: Some NN(j) < 0
353*> -3: NTYPES < 0
354*> -6: THRESH < 0
355*> -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
356*> -14: LDU < 1 or LDU < NMAX.
357*> -26: NWORK too small.
358*> If ZLATMR, CLATMS, or CLATME returns an error code, the
359*> absolute value of it is returned.
360*> If 1, then ZHSEQR could not find all the shifts.
361*> If 2, then the EISPACK code (for small blocks) failed.
362*> If >2, then 30*N iterations were not enough to find an
363*> eigenvalue or to decompose the problem.
364*> Modified.
365*>
366*>-----------------------------------------------------------------------
367*>
368*> Some Local Variables and Parameters:
369*> ---- ----- --------- --- ----------
370*>
371*> ZERO, ONE Real 0 and 1.
372*> MAXTYP The number of types defined.
373*> MTEST The number of tests defined: care must be taken
374*> that (1) the size of RESULT, (2) the number of
375*> tests actually performed, and (3) MTEST agree.
376*> NTEST The number of tests performed on this matrix
377*> so far. This should be less than MTEST, and
378*> equal to it by the last test. It will be less
379*> if any of the routines being tested indicates
380*> that it could not compute the matrices that
381*> would be tested.
382*> NMAX Largest value in NN.
383*> NMATS The number of matrices generated so far.
384*> NERRS The number of tests which have exceeded THRESH
385*> so far (computed by DLAFTS).
386*> COND, CONDS,
387*> IMODE Values to be passed to the matrix generators.
388*> ANORM Norm of A; passed to matrix generators.
389*>
390*> OVFL, UNFL Overflow and underflow thresholds.
391*> ULP, ULPINV Finest relative precision and its inverse.
392*> RTOVFL, RTUNFL,
393*> RTULP, RTULPI Square roots of the previous 4 values.
394*>
395*> The following four arrays decode JTYPE:
396*> KTYPE(j) The general type (1-10) for type "j".
397*> KMODE(j) The MODE value to be passed to the matrix
398*> generator for type "j".
399*> KMAGN(j) The order of magnitude ( O(1),
400*> O(overflow^(1/2) ), O(underflow^(1/2) )
401*> KCONDS(j) Selects whether CONDS is to be 1 or
402*> 1/sqrt(ulp). (0 means irrelevant.)
403*> \endverbatim
404*
405* Authors:
406* ========
407*
408*> \author Univ. of Tennessee
409*> \author Univ. of California Berkeley
410*> \author Univ. of Colorado Denver
411*> \author NAG Ltd.
412*
413*> \ingroup complex16_eig
414*
415* =====================================================================
416 SUBROUTINE zchkhs( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
417 $ NOUNIT, A, LDA, H, T1, T2, U, LDU, Z, UZ, W1,
418 $ W3, EVECTL, EVECTR, EVECTY, EVECTX, UU, TAU,
419 $ WORK, NWORK, RWORK, IWORK, SELECT, RESULT,
420 $ INFO )
421*
422* -- LAPACK test routine --
423* -- LAPACK is a software package provided by Univ. of Tennessee, --
424* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
425*
426* .. Scalar Arguments ..
427 INTEGER INFO, LDA, LDU, NOUNIT, NSIZES, NTYPES, NWORK
428 DOUBLE PRECISION THRESH
429* ..
430* .. Array Arguments ..
431 LOGICAL DOTYPE( * ), SELECT( * )
432 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
433 DOUBLE PRECISION RESULT( 16 ), RWORK( * )
434 COMPLEX*16 A( LDA, * ), EVECTL( LDU, * ),
435 $ evectr( ldu, * ), evectx( ldu, * ),
436 $ evecty( ldu, * ), h( lda, * ), t1( lda, * ),
437 $ t2( lda, * ), tau( * ), u( ldu, * ),
438 $ uu( ldu, * ), uz( ldu, * ), w1( * ), w3( * ),
439 $ work( * ), z( ldu, * )
440* ..
441*
442* =====================================================================
443*
444* .. Parameters ..
445 DOUBLE PRECISION ZERO, ONE
446 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
447 COMPLEX*16 CZERO, CONE
448 PARAMETER ( CZERO = ( 0.0d+0, 0.0d+0 ),
449 $ cone = ( 1.0d+0, 0.0d+0 ) )
450 INTEGER MAXTYP
451 parameter( maxtyp = 21 )
452* ..
453* .. Local Scalars ..
454 LOGICAL BADNN, MATCH
455 INTEGER I, IHI, IINFO, ILO, IMODE, IN, ITYPE, J, JCOL,
456 $ JJ, JSIZE, JTYPE, K, MTYPES, N, N1, NERRS,
457 $ NMATS, NMAX, NTEST, NTESTT
458 DOUBLE PRECISION ANINV, ANORM, COND, CONDS, OVFL, RTOVFL, RTULP,
459 $ rtulpi, rtunfl, temp1, temp2, ulp, ulpinv, unfl
460* ..
461* .. Local Arrays ..
462 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KCONDS( MAXTYP ),
463 $ KMAGN( MAXTYP ), KMODE( MAXTYP ),
464 $ KTYPE( MAXTYP )
465 DOUBLE PRECISION DUMMA( 4 )
466 COMPLEX*16 CDUMMA( 4 )
467* ..
468* .. External Functions ..
469 DOUBLE PRECISION DLAMCH
470 EXTERNAL dlamch
471* ..
472* .. External Subroutines ..
473 EXTERNAL dlafts, dlasum, xerbla, zcopy, zgehrd, zgemm,
476 $ zunghr, zunmhr
477* ..
478* .. Intrinsic Functions ..
479 INTRINSIC abs, dble, max, min, sqrt
480* ..
481* .. Data statements ..
482 DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
483 DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
484 $ 3, 1, 2, 3 /
485 DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
486 $ 1, 5, 5, 5, 4, 3, 1 /
487 DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
488* ..
489* .. Executable Statements ..
490*
491* Check for errors
492*
493 ntestt = 0
494 info = 0
495*
496 badnn = .false.
497 nmax = 0
498 DO 10 j = 1, nsizes
499 nmax = max( nmax, nn( j ) )
500 IF( nn( j ).LT.0 )
501 $ badnn = .true.
502 10 CONTINUE
503*
504* Check for errors
505*
506 IF( nsizes.LT.0 ) THEN
507 info = -1
508 ELSE IF( badnn ) THEN
509 info = -2
510 ELSE IF( ntypes.LT.0 ) THEN
511 info = -3
512 ELSE IF( thresh.LT.zero ) THEN
513 info = -6
514 ELSE IF( lda.LE.1 .OR. lda.LT.nmax ) THEN
515 info = -9
516 ELSE IF( ldu.LE.1 .OR. ldu.LT.nmax ) THEN
517 info = -14
518 ELSE IF( 4*nmax*nmax+2.GT.nwork ) THEN
519 info = -26
520 END IF
521*
522 IF( info.NE.0 ) THEN
523 CALL xerbla( 'ZCHKHS', -info )
524 RETURN
525 END IF
526*
527* Quick return if possible
528*
529 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
530 $ RETURN
531*
532* More important constants
533*
534 unfl = dlamch( 'Safe minimum' )
535 ovfl = dlamch( 'Overflow' )
536 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
537 ulpinv = one / ulp
538 rtunfl = sqrt( unfl )
539 rtovfl = sqrt( ovfl )
540 rtulp = sqrt( ulp )
541 rtulpi = one / rtulp
542*
543* Loop over sizes, types
544*
545 nerrs = 0
546 nmats = 0
547*
548 DO 260 jsize = 1, nsizes
549 n = nn( jsize )
550 IF( n.EQ.0 )
551 $ GO TO 260
552 n1 = max( 1, n )
553 aninv = one / dble( n1 )
554*
555 IF( nsizes.NE.1 ) THEN
556 mtypes = min( maxtyp, ntypes )
557 ELSE
558 mtypes = min( maxtyp+1, ntypes )
559 END IF
560*
561 DO 250 jtype = 1, mtypes
562 IF( .NOT.dotype( jtype ) )
563 $ GO TO 250
564 nmats = nmats + 1
565 ntest = 0
566*
567* Save ISEED in case of an error.
568*
569 DO 20 j = 1, 4
570 ioldsd( j ) = iseed( j )
571 20 CONTINUE
572*
573* Initialize RESULT
574*
575 DO 30 j = 1, 14
576 result( j ) = zero
577 30 CONTINUE
578*
579* Compute "A"
580*
581* Control parameters:
582*
583* KMAGN KCONDS KMODE KTYPE
584* =1 O(1) 1 clustered 1 zero
585* =2 large large clustered 2 identity
586* =3 small exponential Jordan
587* =4 arithmetic diagonal, (w/ eigenvalues)
588* =5 random log hermitian, w/ eigenvalues
589* =6 random general, w/ eigenvalues
590* =7 random diagonal
591* =8 random hermitian
592* =9 random general
593* =10 random triangular
594*
595 IF( mtypes.GT.maxtyp )
596 $ GO TO 100
597*
598 itype = ktype( jtype )
599 imode = kmode( jtype )
600*
601* Compute norm
602*
603 GO TO ( 40, 50, 60 )kmagn( jtype )
604*
605 40 CONTINUE
606 anorm = one
607 GO TO 70
608*
609 50 CONTINUE
610 anorm = ( rtovfl*ulp )*aninv
611 GO TO 70
612*
613 60 CONTINUE
614 anorm = rtunfl*n*ulpinv
615 GO TO 70
616*
617 70 CONTINUE
618*
619 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
620 iinfo = 0
621 cond = ulpinv
622*
623* Special Matrices
624*
625 IF( itype.EQ.1 ) THEN
626*
627* Zero
628*
629 iinfo = 0
630 ELSE IF( itype.EQ.2 ) THEN
631*
632* Identity
633*
634 DO 80 jcol = 1, n
635 a( jcol, jcol ) = anorm
636 80 CONTINUE
637*
638 ELSE IF( itype.EQ.3 ) THEN
639*
640* Jordan Block
641*
642 DO 90 jcol = 1, n
643 a( jcol, jcol ) = anorm
644 IF( jcol.GT.1 )
645 $ a( jcol, jcol-1 ) = one
646 90 CONTINUE
647*
648 ELSE IF( itype.EQ.4 ) THEN
649*
650* Diagonal Matrix, [Eigen]values Specified
651*
652 CALL zlatmr( n, n, 'D', iseed, 'N', work, imode, cond,
653 $ cone, 'T', 'N', work( n+1 ), 1, one,
654 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
655 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
656*
657 ELSE IF( itype.EQ.5 ) THEN
658*
659* Hermitian, eigenvalues specified
660*
661 CALL zlatms( n, n, 'D', iseed, 'H', rwork, imode, cond,
662 $ anorm, n, n, 'N', a, lda, work, iinfo )
663*
664 ELSE IF( itype.EQ.6 ) THEN
665*
666* General, eigenvalues specified
667*
668 IF( kconds( jtype ).EQ.1 ) THEN
669 conds = one
670 ELSE IF( kconds( jtype ).EQ.2 ) THEN
671 conds = rtulpi
672 ELSE
673 conds = zero
674 END IF
675*
676 CALL zlatme( n, 'D', iseed, work, imode, cond, cone,
677 $ 'T', 'T', 'T', rwork, 4, conds, n, n, anorm,
678 $ a, lda, work( n+1 ), iinfo )
679*
680 ELSE IF( itype.EQ.7 ) THEN
681*
682* Diagonal, random eigenvalues
683*
684 CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
685 $ 'T', 'N', work( n+1 ), 1, one,
686 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
687 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
688*
689 ELSE IF( itype.EQ.8 ) THEN
690*
691* Hermitian, random eigenvalues
692*
693 CALL zlatmr( n, n, 'D', iseed, 'H', work, 6, one, cone,
694 $ 'T', 'N', work( n+1 ), 1, one,
695 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
696 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
697*
698 ELSE IF( itype.EQ.9 ) THEN
699*
700* General, random eigenvalues
701*
702 CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
703 $ 'T', 'N', work( n+1 ), 1, one,
704 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
705 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
706*
707 ELSE IF( itype.EQ.10 ) THEN
708*
709* Triangular, random eigenvalues
710*
711 CALL zlatmr( n, n, 'D', iseed, 'N', work, 6, one, cone,
712 $ 'T', 'N', work( n+1 ), 1, one,
713 $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
714 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
715*
716 ELSE
717*
718 iinfo = 1
719 END IF
720*
721 IF( iinfo.NE.0 ) THEN
722 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
723 $ ioldsd
724 info = abs( iinfo )
725 RETURN
726 END IF
727*
728 100 CONTINUE
729*
730* Call ZGEHRD to compute H and U, do tests.
731*
732 CALL zlacpy( ' ', n, n, a, lda, h, lda )
733 ntest = 1
734*
735 ilo = 1
736 ihi = n
737*
738 CALL zgehrd( n, ilo, ihi, h, lda, work, work( n+1 ),
739 $ nwork-n, iinfo )
740*
741 IF( iinfo.NE.0 ) THEN
742 result( 1 ) = ulpinv
743 WRITE( nounit, fmt = 9999 )'ZGEHRD', iinfo, n, jtype,
744 $ ioldsd
745 info = abs( iinfo )
746 GO TO 240
747 END IF
748*
749 DO 120 j = 1, n - 1
750 uu( j+1, j ) = czero
751 DO 110 i = j + 2, n
752 u( i, j ) = h( i, j )
753 uu( i, j ) = h( i, j )
754 h( i, j ) = czero
755 110 CONTINUE
756 120 CONTINUE
757 CALL zcopy( n-1, work, 1, tau, 1 )
758 CALL zunghr( n, ilo, ihi, u, ldu, work, work( n+1 ),
759 $ nwork-n, iinfo )
760 ntest = 2
761*
762 CALL zhst01( n, ilo, ihi, a, lda, h, lda, u, ldu, work,
763 $ nwork, rwork, result( 1 ) )
764*
765* Call ZHSEQR to compute T1, T2 and Z, do tests.
766*
767* Eigenvalues only (W3)
768*
769 CALL zlacpy( ' ', n, n, h, lda, t2, lda )
770 ntest = 3
771 result( 3 ) = ulpinv
772*
773 CALL zhseqr( 'E', 'N', n, ilo, ihi, t2, lda, w3, uz, ldu,
774 $ work, nwork, iinfo )
775 IF( iinfo.NE.0 ) THEN
776 WRITE( nounit, fmt = 9999 )'ZHSEQR(E)', iinfo, n, jtype,
777 $ ioldsd
778 IF( iinfo.LE.n+2 ) THEN
779 info = abs( iinfo )
780 GO TO 240
781 END IF
782 END IF
783*
784* Eigenvalues (W1) and Full Schur Form (T2)
785*
786 CALL zlacpy( ' ', n, n, h, lda, t2, lda )
787*
788 CALL zhseqr( 'S', 'N', n, ilo, ihi, t2, lda, w1, uz, ldu,
789 $ work, nwork, iinfo )
790 IF( iinfo.NE.0 .AND. iinfo.LE.n+2 ) THEN
791 WRITE( nounit, fmt = 9999 )'ZHSEQR(S)', iinfo, n, jtype,
792 $ ioldsd
793 info = abs( iinfo )
794 GO TO 240
795 END IF
796*
797* Eigenvalues (W1), Schur Form (T1), and Schur Vectors (UZ)
798*
799 CALL zlacpy( ' ', n, n, h, lda, t1, lda )
800 CALL zlacpy( ' ', n, n, u, ldu, uz, ldu )
801*
802 CALL zhseqr( 'S', 'V', n, ilo, ihi, t1, lda, w1, uz, ldu,
803 $ work, nwork, iinfo )
804 IF( iinfo.NE.0 .AND. iinfo.LE.n+2 ) THEN
805 WRITE( nounit, fmt = 9999 )'ZHSEQR(V)', iinfo, n, jtype,
806 $ ioldsd
807 info = abs( iinfo )
808 GO TO 240
809 END IF
810*
811* Compute Z = U' UZ
812*
813 CALL zgemm( 'C', 'N', n, n, n, cone, u, ldu, uz, ldu, czero,
814 $ z, ldu )
815 ntest = 8
816*
817* Do Tests 3: | H - Z T Z' | / ( |H| n ulp )
818* and 4: | I - Z Z' | / ( n ulp )
819*
820 CALL zhst01( n, ilo, ihi, h, lda, t1, lda, z, ldu, work,
821 $ nwork, rwork, result( 3 ) )
822*
823* Do Tests 5: | A - UZ T (UZ)' | / ( |A| n ulp )
824* and 6: | I - UZ (UZ)' | / ( n ulp )
825*
826 CALL zhst01( n, ilo, ihi, a, lda, t1, lda, uz, ldu, work,
827 $ nwork, rwork, result( 5 ) )
828*
829* Do Test 7: | T2 - T1 | / ( |T| n ulp )
830*
831 CALL zget10( n, n, t2, lda, t1, lda, work, rwork,
832 $ result( 7 ) )
833*
834* Do Test 8: | W3 - W1 | / ( max(|W1|,|W3|) ulp )
835*
836 temp1 = zero
837 temp2 = zero
838 DO 130 j = 1, n
839 temp1 = max( temp1, abs( w1( j ) ), abs( w3( j ) ) )
840 temp2 = max( temp2, abs( w1( j )-w3( j ) ) )
841 130 CONTINUE
842*
843 result( 8 ) = temp2 / max( unfl, ulp*max( temp1, temp2 ) )
844*
845* Compute the Left and Right Eigenvectors of T
846*
847* Compute the Right eigenvector Matrix:
848*
849 ntest = 9
850 result( 9 ) = ulpinv
851*
852* Select every other eigenvector
853*
854 DO 140 j = 1, n
855 SELECT( j ) = .false.
856 140 CONTINUE
857 DO 150 j = 1, n, 2
858 SELECT( j ) = .true.
859 150 CONTINUE
860 CALL ztrevc( 'Right', 'All', SELECT, n, t1, lda, cdumma,
861 $ ldu, evectr, ldu, n, in, work, rwork, iinfo )
862 IF( iinfo.NE.0 ) THEN
863 WRITE( nounit, fmt = 9999 )'ZTREVC(R,A)', iinfo, n,
864 $ jtype, ioldsd
865 info = abs( iinfo )
866 GO TO 240
867 END IF
868*
869* Test 9: | TR - RW | / ( |T| |R| ulp )
870*
871 CALL zget22( 'N', 'N', 'N', n, t1, lda, evectr, ldu, w1,
872 $ work, rwork, dumma( 1 ) )
873 result( 9 ) = dumma( 1 )
874 IF( dumma( 2 ).GT.thresh ) THEN
875 WRITE( nounit, fmt = 9998 )'Right', 'ZTREVC',
876 $ dumma( 2 ), n, jtype, ioldsd
877 END IF
878*
879* Compute selected right eigenvectors and confirm that
880* they agree with previous right eigenvectors
881*
882 CALL ztrevc( 'Right', 'Some', SELECT, n, t1, lda, cdumma,
883 $ ldu, evectl, ldu, n, in, work, rwork, iinfo )
884 IF( iinfo.NE.0 ) THEN
885 WRITE( nounit, fmt = 9999 )'ZTREVC(R,S)', iinfo, n,
886 $ jtype, ioldsd
887 info = abs( iinfo )
888 GO TO 240
889 END IF
890*
891 k = 1
892 match = .true.
893 DO 170 j = 1, n
894 IF( SELECT( j ) ) THEN
895 DO 160 jj = 1, n
896 IF( evectr( jj, j ).NE.evectl( jj, k ) ) THEN
897 match = .false.
898 GO TO 180
899 END IF
900 160 CONTINUE
901 k = k + 1
902 END IF
903 170 CONTINUE
904 180 CONTINUE
905 IF( .NOT.match )
906 $ WRITE( nounit, fmt = 9997 )'Right', 'ZTREVC', n, jtype,
907 $ ioldsd
908*
909* Compute the Left eigenvector Matrix:
910*
911 ntest = 10
912 result( 10 ) = ulpinv
913 CALL ztrevc( 'Left', 'All', SELECT, n, t1, lda, evectl, ldu,
914 $ cdumma, ldu, n, in, work, rwork, iinfo )
915 IF( iinfo.NE.0 ) THEN
916 WRITE( nounit, fmt = 9999 )'ZTREVC(L,A)', iinfo, n,
917 $ jtype, ioldsd
918 info = abs( iinfo )
919 GO TO 240
920 END IF
921*
922* Test 10: | LT - WL | / ( |T| |L| ulp )
923*
924 CALL zget22( 'C', 'N', 'C', n, t1, lda, evectl, ldu, w1,
925 $ work, rwork, dumma( 3 ) )
926 result( 10 ) = dumma( 3 )
927 IF( dumma( 4 ).GT.thresh ) THEN
928 WRITE( nounit, fmt = 9998 )'Left', 'ZTREVC', dumma( 4 ),
929 $ n, jtype, ioldsd
930 END IF
931*
932* Compute selected left eigenvectors and confirm that
933* they agree with previous left eigenvectors
934*
935 CALL ztrevc( 'Left', 'Some', SELECT, n, t1, lda, evectr,
936 $ ldu, cdumma, ldu, n, in, work, rwork, iinfo )
937 IF( iinfo.NE.0 ) THEN
938 WRITE( nounit, fmt = 9999 )'ZTREVC(L,S)', iinfo, n,
939 $ jtype, ioldsd
940 info = abs( iinfo )
941 GO TO 240
942 END IF
943*
944 k = 1
945 match = .true.
946 DO 200 j = 1, n
947 IF( SELECT( j ) ) THEN
948 DO 190 jj = 1, n
949 IF( evectl( jj, j ).NE.evectr( jj, k ) ) THEN
950 match = .false.
951 GO TO 210
952 END IF
953 190 CONTINUE
954 k = k + 1
955 END IF
956 200 CONTINUE
957 210 CONTINUE
958 IF( .NOT.match )
959 $ WRITE( nounit, fmt = 9997 )'Left', 'ZTREVC', n, jtype,
960 $ ioldsd
961*
962* Call ZHSEIN for Right eigenvectors of H, do test 11
963*
964 ntest = 11
965 result( 11 ) = ulpinv
966 DO 220 j = 1, n
967 SELECT( j ) = .true.
968 220 CONTINUE
969*
970 CALL zhsein( 'Right', 'Qr', 'Ninitv', SELECT, n, h, lda, w3,
971 $ cdumma, ldu, evectx, ldu, n1, in, work, rwork,
972 $ iwork, iwork, iinfo )
973 IF( iinfo.NE.0 ) THEN
974 WRITE( nounit, fmt = 9999 )'ZHSEIN(R)', iinfo, n, jtype,
975 $ ioldsd
976 info = abs( iinfo )
977 IF( iinfo.LT.0 )
978 $ GO TO 240
979 ELSE
980*
981* Test 11: | HX - XW | / ( |H| |X| ulp )
982*
983* (from inverse iteration)
984*
985 CALL zget22( 'N', 'N', 'N', n, h, lda, evectx, ldu, w3,
986 $ work, rwork, dumma( 1 ) )
987 IF( dumma( 1 ).LT.ulpinv )
988 $ result( 11 ) = dumma( 1 )*aninv
989 IF( dumma( 2 ).GT.thresh ) THEN
990 WRITE( nounit, fmt = 9998 )'Right', 'ZHSEIN',
991 $ dumma( 2 ), n, jtype, ioldsd
992 END IF
993 END IF
994*
995* Call ZHSEIN for Left eigenvectors of H, do test 12
996*
997 ntest = 12
998 result( 12 ) = ulpinv
999 DO 230 j = 1, n
1000 SELECT( j ) = .true.
1001 230 CONTINUE
1002*
1003 CALL zhsein( 'Left', 'Qr', 'Ninitv', SELECT, n, h, lda, w3,
1004 $ evecty, ldu, cdumma, ldu, n1, in, work, rwork,
1005 $ iwork, iwork, iinfo )
1006 IF( iinfo.NE.0 ) THEN
1007 WRITE( nounit, fmt = 9999 )'ZHSEIN(L)', iinfo, n, jtype,
1008 $ ioldsd
1009 info = abs( iinfo )
1010 IF( iinfo.LT.0 )
1011 $ GO TO 240
1012 ELSE
1013*
1014* Test 12: | YH - WY | / ( |H| |Y| ulp )
1015*
1016* (from inverse iteration)
1017*
1018 CALL zget22( 'C', 'N', 'C', n, h, lda, evecty, ldu, w3,
1019 $ work, rwork, dumma( 3 ) )
1020 IF( dumma( 3 ).LT.ulpinv )
1021 $ result( 12 ) = dumma( 3 )*aninv
1022 IF( dumma( 4 ).GT.thresh ) THEN
1023 WRITE( nounit, fmt = 9998 )'Left', 'ZHSEIN',
1024 $ dumma( 4 ), n, jtype, ioldsd
1025 END IF
1026 END IF
1027*
1028* Call ZUNMHR for Right eigenvectors of A, do test 13
1029*
1030 ntest = 13
1031 result( 13 ) = ulpinv
1032*
1033 CALL zunmhr( 'Left', 'No transpose', n, n, ilo, ihi, uu,
1034 $ ldu, tau, evectx, ldu, work, nwork, iinfo )
1035 IF( iinfo.NE.0 ) THEN
1036 WRITE( nounit, fmt = 9999 )'ZUNMHR(L)', iinfo, n, jtype,
1037 $ ioldsd
1038 info = abs( iinfo )
1039 IF( iinfo.LT.0 )
1040 $ GO TO 240
1041 ELSE
1042*
1043* Test 13: | AX - XW | / ( |A| |X| ulp )
1044*
1045* (from inverse iteration)
1046*
1047 CALL zget22( 'N', 'N', 'N', n, a, lda, evectx, ldu, w3,
1048 $ work, rwork, dumma( 1 ) )
1049 IF( dumma( 1 ).LT.ulpinv )
1050 $ result( 13 ) = dumma( 1 )*aninv
1051 END IF
1052*
1053* Call ZUNMHR for Left eigenvectors of A, do test 14
1054*
1055 ntest = 14
1056 result( 14 ) = ulpinv
1057*
1058 CALL zunmhr( 'Left', 'No transpose', n, n, ilo, ihi, uu,
1059 $ ldu, tau, evecty, ldu, work, nwork, iinfo )
1060 IF( iinfo.NE.0 ) THEN
1061 WRITE( nounit, fmt = 9999 )'ZUNMHR(L)', iinfo, n, jtype,
1062 $ ioldsd
1063 info = abs( iinfo )
1064 IF( iinfo.LT.0 )
1065 $ GO TO 240
1066 ELSE
1067*
1068* Test 14: | YA - WY | / ( |A| |Y| ulp )
1069*
1070* (from inverse iteration)
1071*
1072 CALL zget22( 'C', 'N', 'C', n, a, lda, evecty, ldu, w3,
1073 $ work, rwork, dumma( 3 ) )
1074 IF( dumma( 3 ).LT.ulpinv )
1075 $ result( 14 ) = dumma( 3 )*aninv
1076 END IF
1077*
1078* Compute Left and Right Eigenvectors of A
1079*
1080* Compute a Right eigenvector matrix:
1081*
1082 ntest = 15
1083 result( 15 ) = ulpinv
1084*
1085 CALL zlacpy( ' ', n, n, uz, ldu, evectr, ldu )
1086*
1087 CALL ztrevc3( 'Right', 'Back', SELECT, n, t1, lda, cdumma,
1088 $ ldu, evectr, ldu, n, in, work, nwork, rwork,
1089 $ n, iinfo )
1090 IF( iinfo.NE.0 ) THEN
1091 WRITE( nounit, fmt = 9999 )'ZTREVC3(R,B)', iinfo, n,
1092 $ jtype, ioldsd
1093 info = abs( iinfo )
1094 GO TO 250
1095 END IF
1096*
1097* Test 15: | AR - RW | / ( |A| |R| ulp )
1098*
1099* (from Schur decomposition)
1100*
1101 CALL zget22( 'N', 'N', 'N', n, a, lda, evectr, ldu, w1,
1102 $ work, rwork, dumma( 1 ) )
1103 result( 15 ) = dumma( 1 )
1104 IF( dumma( 2 ).GT.thresh ) THEN
1105 WRITE( nounit, fmt = 9998 )'Right', 'ZTREVC3',
1106 $ dumma( 2 ), n, jtype, ioldsd
1107 END IF
1108*
1109* Compute a Left eigenvector matrix:
1110*
1111 ntest = 16
1112 result( 16 ) = ulpinv
1113*
1114 CALL zlacpy( ' ', n, n, uz, ldu, evectl, ldu )
1115*
1116 CALL ztrevc3( 'Left', 'Back', SELECT, n, t1, lda, evectl,
1117 $ ldu, cdumma, ldu, n, in, work, nwork, rwork,
1118 $ n, iinfo )
1119 IF( iinfo.NE.0 ) THEN
1120 WRITE( nounit, fmt = 9999 )'ZTREVC3(L,B)', iinfo, n,
1121 $ jtype, ioldsd
1122 info = abs( iinfo )
1123 GO TO 250
1124 END IF
1125*
1126* Test 16: | LA - WL | / ( |A| |L| ulp )
1127*
1128* (from Schur decomposition)
1129*
1130 CALL zget22( 'Conj', 'N', 'Conj', n, a, lda, evectl, ldu,
1131 $ w1, work, rwork, dumma( 3 ) )
1132 result( 16 ) = dumma( 3 )
1133 IF( dumma( 4 ).GT.thresh ) THEN
1134 WRITE( nounit, fmt = 9998 )'Left', 'ZTREVC3', dumma( 4 ),
1135 $ n, jtype, ioldsd
1136 END IF
1137*
1138* End of Loop -- Check for RESULT(j) > THRESH
1139*
1140 240 CONTINUE
1141*
1142 ntestt = ntestt + ntest
1143 CALL dlafts( 'ZHS', n, n, jtype, ntest, result, ioldsd,
1144 $ thresh, nounit, nerrs )
1145*
1146 250 CONTINUE
1147 260 CONTINUE
1148*
1149* Summary
1150*
1151 CALL dlasum( 'ZHS', nounit, nerrs, ntestt )
1152*
1153 RETURN
1154*
1155 9999 FORMAT( ' ZCHKHS: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1156 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1157 9998 FORMAT( ' ZCHKHS: ', a, ' Eigenvectors from ', a, ' incorrectly ',
1158 $ 'normalized.', / ' Bits of error=', 0p, g10.3, ',', 9x,
1159 $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
1160 $ ')' )
1161 9997 FORMAT( ' ZCHKHS: Selected ', a, ' Eigenvectors from ', a,
1162 $ ' do not match other eigenvectors ', 9x, 'N=', i6,
1163 $ ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1164*
1165* End of ZCHKHS
1166*
1167 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlafts(type, m, n, imat, ntests, result, iseed, thresh, iounit, ie)
DLAFTS
Definition dlafts.f:99
subroutine dlasum(type, iounit, ie, nrun)
DLASUM
Definition dlasum.f:43
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
Definition zgehrd.f:167
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zhsein(side, eigsrc, initv, select, n, h, ldh, w, vl, ldvl, vr, ldvr, mm, m, work, rwork, ifaill, ifailr, info)
ZHSEIN
Definition zhsein.f:245
subroutine zhseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
ZHSEQR
Definition zhseqr.f:299
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine ztrevc3(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, lwork, rwork, lrwork, info)
ZTREVC3
Definition ztrevc3.f:244
subroutine ztrevc(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
ZTREVC
Definition ztrevc.f:218
subroutine zunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZUNGHR
Definition zunghr.f:126
subroutine zunmhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
ZUNMHR
Definition zunmhr.f:178
subroutine zchkhs(nsizes, nn, ntypes, dotype, iseed, thresh, nounit, a, lda, h, t1, t2, u, ldu, z, uz, w1, w3, evectl, evectr, evecty, evectx, uu, tau, work, nwork, rwork, iwork, select, result, info)
ZCHKHS
Definition zchkhs.f:421
subroutine zget10(m, n, a, lda, b, ldb, work, rwork, result)
ZGET10
Definition zget10.f:99
subroutine zget22(transa, transe, transw, n, a, lda, e, lde, w, work, rwork, result)
ZGET22
Definition zget22.f:144
subroutine zhst01(n, ilo, ihi, a, lda, h, ldh, q, ldq, work, lwork, rwork, result)
ZHST01
Definition zhst01.f:140
subroutine zlatme(n, dist, iseed, d, mode, cond, dmax, rsign, upper, sim, ds, modes, conds, kl, ku, anorm, a, lda, work, info)
ZLATME
Definition zlatme.f:301
subroutine zlatmr(m, n, dist, iseed, sym, d, mode, cond, dmax, rsign, grade, dl, model, condl, dr, moder, condr, pivtng, ipivot, kl, ku, sparse, anorm, pack, a, lda, iwork, info)
ZLATMR
Definition zlatmr.f:490
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332