LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvst2stg.f
Go to the documentation of this file.
1*> \brief \b ZDRVST2STG
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 ZDRVST2STG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12* NOUNIT, A, LDA, D1, D2, D3, WA1, WA2, WA3, U,
13* LDU, V, TAU, Z, WORK, LWORK, RWORK, LRWORK,
14* IWORK, LIWORK, RESULT, INFO )
15*
16* .. Scalar Arguments ..
17* INTEGER INFO, LDA, LDU, LIWORK, LRWORK, LWORK, NOUNIT,
18* $ NSIZES, NTYPES
19* DOUBLE PRECISION THRESH
20* ..
21* .. Array Arguments ..
22* LOGICAL DOTYPE( * )
23* INTEGER ISEED( 4 ), IWORK( * ), NN( * )
24* DOUBLE PRECISION D1( * ), D2( * ), D3( * ), RESULT( * ),
25* $ RWORK( * ), WA1( * ), WA2( * ), WA3( * )
26* COMPLEX*16 A( LDA, * ), TAU( * ), U( LDU, * ),
27* $ V( LDU, * ), WORK( * ), Z( LDU, * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> ZDRVST2STG checks the Hermitian eigenvalue problem drivers.
37*>
38*> ZHEEVD computes all eigenvalues and, optionally,
39*> eigenvectors of a complex Hermitian matrix,
40*> using a divide-and-conquer algorithm.
41*>
42*> ZHEEVX computes selected eigenvalues and, optionally,
43*> eigenvectors of a complex Hermitian matrix.
44*>
45*> ZHEEVR computes selected eigenvalues and, optionally,
46*> eigenvectors of a complex Hermitian matrix
47*> using the Relatively Robust Representation where it can.
48*>
49*> ZHPEVD computes all eigenvalues and, optionally,
50*> eigenvectors of a complex Hermitian matrix in packed
51*> storage, using a divide-and-conquer algorithm.
52*>
53*> ZHPEVX computes selected eigenvalues and, optionally,
54*> eigenvectors of a complex Hermitian matrix in packed
55*> storage.
56*>
57*> ZHBEVD computes all eigenvalues and, optionally,
58*> eigenvectors of a complex Hermitian band matrix,
59*> using a divide-and-conquer algorithm.
60*>
61*> ZHBEVX computes selected eigenvalues and, optionally,
62*> eigenvectors of a complex Hermitian band matrix.
63*>
64*> ZHEEV computes all eigenvalues and, optionally,
65*> eigenvectors of a complex Hermitian matrix.
66*>
67*> ZHPEV computes all eigenvalues and, optionally,
68*> eigenvectors of a complex Hermitian matrix in packed
69*> storage.
70*>
71*> ZHBEV computes all eigenvalues and, optionally,
72*> eigenvectors of a complex Hermitian band matrix.
73*>
74*> When ZDRVST2STG is called, a number of matrix "sizes" ("n's") and a
75*> number of matrix "types" are specified. For each size ("n")
76*> and each type of matrix, one matrix will be generated and used
77*> to test the appropriate drivers. For each matrix and each
78*> driver routine called, the following tests will be performed:
79*>
80*> (1) | A - Z D Z' | / ( |A| n ulp )
81*>
82*> (2) | I - Z Z' | / ( n ulp )
83*>
84*> (3) | D1 - D2 | / ( |D1| ulp )
85*>
86*> where Z is the matrix of eigenvectors returned when the
87*> eigenvector option is given and D1 and D2 are the eigenvalues
88*> returned with and without the eigenvector option.
89*>
90*> The "sizes" are specified by an array NN(1:NSIZES); the value of
91*> each element NN(j) specifies one size.
92*> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
93*> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
94*> Currently, the list of possible types is:
95*>
96*> (1) The zero matrix.
97*> (2) The identity matrix.
98*>
99*> (3) A diagonal matrix with evenly spaced entries
100*> 1, ..., ULP and random signs.
101*> (ULP = (first number larger than 1) - 1 )
102*> (4) A diagonal matrix with geometrically spaced entries
103*> 1, ..., ULP and random signs.
104*> (5) A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
105*> and random signs.
106*>
107*> (6) Same as (4), but multiplied by SQRT( overflow threshold )
108*> (7) Same as (4), but multiplied by SQRT( underflow threshold )
109*>
110*> (8) A matrix of the form U* D U, where U is unitary and
111*> D has evenly spaced entries 1, ..., ULP with random signs
112*> on the diagonal.
113*>
114*> (9) A matrix of the form U* D U, where U is unitary and
115*> D has geometrically spaced entries 1, ..., ULP with random
116*> signs on the diagonal.
117*>
118*> (10) A matrix of the form U* D U, where U is unitary and
119*> D has "clustered" entries 1, ULP,..., ULP with random
120*> signs on the diagonal.
121*>
122*> (11) Same as (8), but multiplied by SQRT( overflow threshold )
123*> (12) Same as (8), but multiplied by SQRT( underflow threshold )
124*>
125*> (13) Symmetric matrix with random entries chosen from (-1,1).
126*> (14) Same as (13), but multiplied by SQRT( overflow threshold )
127*> (15) Same as (13), but multiplied by SQRT( underflow threshold )
128*> (16) A band matrix with half bandwidth randomly chosen between
129*> 0 and N-1, with evenly spaced eigenvalues 1, ..., ULP
130*> with random signs.
131*> (17) Same as (16), but multiplied by SQRT( overflow threshold )
132*> (18) Same as (16), but multiplied by SQRT( underflow threshold )
133*> \endverbatim
134*
135* Arguments:
136* ==========
137*
138*> \verbatim
139*> NSIZES INTEGER
140*> The number of sizes of matrices to use. If it is zero,
141*> ZDRVST2STG does nothing. It must be at least zero.
142*> Not modified.
143*>
144*> NN INTEGER array, dimension (NSIZES)
145*> An array containing the sizes to be used for the matrices.
146*> Zero values will be skipped. The values must be at least
147*> zero.
148*> Not modified.
149*>
150*> NTYPES INTEGER
151*> The number of elements in DOTYPE. If it is zero, ZDRVST2STG
152*> does nothing. It must be at least zero. If it is MAXTYP+1
153*> and NSIZES is 1, then an additional type, MAXTYP+1 is
154*> defined, which is to use whatever matrix is in A. This
155*> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
156*> DOTYPE(MAXTYP+1) is .TRUE. .
157*> Not modified.
158*>
159*> DOTYPE LOGICAL array, dimension (NTYPES)
160*> If DOTYPE(j) is .TRUE., then for each size in NN a
161*> matrix of that size and of type j will be generated.
162*> If NTYPES is smaller than the maximum number of types
163*> defined (PARAMETER MAXTYP), then types NTYPES+1 through
164*> MAXTYP will not be generated. If NTYPES is larger
165*> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
166*> will be ignored.
167*> Not modified.
168*>
169*> ISEED INTEGER array, dimension (4)
170*> On entry ISEED specifies the seed of the random number
171*> generator. The array elements should be between 0 and 4095;
172*> if not they will be reduced mod 4096. Also, ISEED(4) must
173*> be odd. The random number generator uses a linear
174*> congruential sequence limited to small integers, and so
175*> should produce machine independent random numbers. The
176*> values of ISEED are changed on exit, and can be used in the
177*> next call to ZDRVST2STG to continue the same random number
178*> sequence.
179*> Modified.
180*>
181*> THRESH DOUBLE PRECISION
182*> A test will count as "failed" if the "error", computed as
183*> described above, exceeds THRESH. Note that the error
184*> is scaled to be O(1), so THRESH should be a reasonably
185*> small multiple of 1, e.g., 10 or 100. In particular,
186*> it should not depend on the precision (single vs. double)
187*> or the size of the matrix. It must be at least zero.
188*> Not modified.
189*>
190*> NOUNIT INTEGER
191*> The FORTRAN unit number for printing out error messages
192*> (e.g., if a routine returns IINFO not equal to 0.)
193*> Not modified.
194*>
195*> A COMPLEX*16 array, dimension (LDA , max(NN))
196*> Used to hold the matrix whose eigenvalues are to be
197*> computed. On exit, A contains the last matrix actually
198*> used.
199*> Modified.
200*>
201*> LDA INTEGER
202*> The leading dimension of A. It must be at
203*> least 1 and at least max( NN ).
204*> Not modified.
205*>
206*> D1 DOUBLE PRECISION array, dimension (max(NN))
207*> The eigenvalues of A, as computed by ZSTEQR simultaneously
208*> with Z. On exit, the eigenvalues in D1 correspond with the
209*> matrix in A.
210*> Modified.
211*>
212*> D2 DOUBLE PRECISION array, dimension (max(NN))
213*> The eigenvalues of A, as computed by ZSTEQR if Z is not
214*> computed. On exit, the eigenvalues in D2 correspond with
215*> the matrix in A.
216*> Modified.
217*>
218*> D3 DOUBLE PRECISION array, dimension (max(NN))
219*> The eigenvalues of A, as computed by DSTERF. On exit, the
220*> eigenvalues in D3 correspond with the matrix in A.
221*> Modified.
222*>
223*> WA1 DOUBLE PRECISION array, dimension
224*>
225*> WA2 DOUBLE PRECISION array, dimension
226*>
227*> WA3 DOUBLE PRECISION array, dimension
228*>
229*> U COMPLEX*16 array, dimension (LDU, max(NN))
230*> The unitary matrix computed by ZHETRD + ZUNGC3.
231*> Modified.
232*>
233*> LDU INTEGER
234*> The leading dimension of U, Z, and V. It must be at
235*> least 1 and at least max( NN ).
236*> Not modified.
237*>
238*> V COMPLEX*16 array, dimension (LDU, max(NN))
239*> The Housholder vectors computed by ZHETRD in reducing A to
240*> tridiagonal form.
241*> Modified.
242*>
243*> TAU COMPLEX*16 array, dimension (max(NN))
244*> The Householder factors computed by ZHETRD in reducing A
245*> to tridiagonal form.
246*> Modified.
247*>
248*> Z COMPLEX*16 array, dimension (LDU, max(NN))
249*> The unitary matrix of eigenvectors computed by ZHEEVD,
250*> ZHEEVX, ZHPEVD, CHPEVX, ZHBEVD, and CHBEVX.
251*> Modified.
252*>
253*> WORK - COMPLEX*16 array of dimension ( LWORK )
254*> Workspace.
255*> Modified.
256*>
257*> LWORK - INTEGER
258*> The number of entries in WORK. This must be at least
259*> 2*max( NN(j), 2 )**2.
260*> Not modified.
261*>
262*> RWORK DOUBLE PRECISION array, dimension (3*max(NN))
263*> Workspace.
264*> Modified.
265*>
266*> LRWORK - INTEGER
267*> The number of entries in RWORK.
268*>
269*> IWORK INTEGER array, dimension (6*max(NN))
270*> Workspace.
271*> Modified.
272*>
273*> LIWORK - INTEGER
274*> The number of entries in IWORK.
275*>
276*> RESULT DOUBLE PRECISION array, dimension (??)
277*> The values computed by the tests described above.
278*> The values are currently limited to 1/ulp, to avoid
279*> overflow.
280*> Modified.
281*>
282*> INFO INTEGER
283*> If 0, then everything ran OK.
284*> -1: NSIZES < 0
285*> -2: Some NN(j) < 0
286*> -3: NTYPES < 0
287*> -5: THRESH < 0
288*> -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
289*> -16: LDU < 1 or LDU < NMAX.
290*> -21: LWORK too small.
291*> If DLATMR, SLATMS, ZHETRD, DORGC3, ZSTEQR, DSTERF,
292*> or DORMC2 returns an error code, the
293*> absolute value of it is returned.
294*> Modified.
295*>
296*>-----------------------------------------------------------------------
297*>
298*> Some Local Variables and Parameters:
299*> ---- ----- --------- --- ----------
300*> ZERO, ONE Real 0 and 1.
301*> MAXTYP The number of types defined.
302*> NTEST The number of tests performed, or which can
303*> be performed so far, for the current matrix.
304*> NTESTT The total number of tests performed so far.
305*> NMAX Largest value in NN.
306*> NMATS The number of matrices generated so far.
307*> NERRS The number of tests which have exceeded THRESH
308*> so far (computed by DLAFTS).
309*> COND, IMODE Values to be passed to the matrix generators.
310*> ANORM Norm of A; passed to matrix generators.
311*>
312*> OVFL, UNFL Overflow and underflow thresholds.
313*> ULP, ULPINV Finest relative precision and its inverse.
314*> RTOVFL, RTUNFL Square roots of the previous 2 values.
315*> The following four arrays decode JTYPE:
316*> KTYPE(j) The general type (1-10) for type "j".
317*> KMODE(j) The MODE value to be passed to the matrix
318*> generator for type "j".
319*> KMAGN(j) The order of magnitude ( O(1),
320*> O(overflow^(1/2) ), O(underflow^(1/2) )
321*> \endverbatim
322*
323* Authors:
324* ========
325*
326*> \author Univ. of Tennessee
327*> \author Univ. of California Berkeley
328*> \author Univ. of Colorado Denver
329*> \author NAG Ltd.
330*
331*> \ingroup complex16_eig
332*
333* =====================================================================
334 SUBROUTINE zdrvst2stg( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
335 $ NOUNIT, A, LDA, D1, D2, D3, WA1, WA2, WA3, U,
336 $ LDU, V, TAU, Z, WORK, LWORK, RWORK, LRWORK,
337 $ IWORK, LIWORK, RESULT, INFO )
338*
339* -- LAPACK test routine --
340* -- LAPACK is a software package provided by Univ. of Tennessee, --
341* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
342*
343* .. Scalar Arguments ..
344 INTEGER INFO, LDA, LDU, LIWORK, LRWORK, LWORK, NOUNIT,
345 $ NSIZES, NTYPES
346 DOUBLE PRECISION THRESH
347* ..
348* .. Array Arguments ..
349 LOGICAL DOTYPE( * )
350 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
351 DOUBLE PRECISION D1( * ), D2( * ), D3( * ), RESULT( * ),
352 $ rwork( * ), wa1( * ), wa2( * ), wa3( * )
353 COMPLEX*16 A( LDA, * ), TAU( * ), U( LDU, * ),
354 $ v( ldu, * ), work( * ), z( ldu, * )
355* ..
356*
357* =====================================================================
358*
359*
360* .. Parameters ..
361 DOUBLE PRECISION ZERO, ONE, TWO, TEN
362 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0, two = 2.0d+0,
363 $ ten = 10.0d+0 )
364 DOUBLE PRECISION HALF
365 parameter( half = one / two )
366 COMPLEX*16 CZERO, CONE
367 parameter( czero = ( 0.0d+0, 0.0d+0 ),
368 $ cone = ( 1.0d+0, 0.0d+0 ) )
369 INTEGER MAXTYP
370 parameter( maxtyp = 18 )
371* ..
372* .. Local Scalars ..
373 LOGICAL BADNN
374 CHARACTER UPLO
375 INTEGER I, IDIAG, IHBW, IINFO, IL, IMODE, INDWRK, INDX,
376 $ irow, itemp, itype, iu, iuplo, j, j1, j2, jcol,
377 $ jsize, jtype, kd, lgn, liwedc, lrwedc, lwedc,
378 $ m, m2, m3, mtypes, n, nerrs, nmats, nmax,
379 $ ntest, ntestt
380 DOUBLE PRECISION ABSTOL, ANINV, ANORM, COND, OVFL, RTOVFL,
381 $ RTUNFL, TEMP1, TEMP2, TEMP3, ULP, ULPINV, UNFL,
382 $ VL, VU
383* ..
384* .. Local Arrays ..
385 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), ISEED2( 4 ),
386 $ ISEED3( 4 ), KMAGN( MAXTYP ), KMODE( MAXTYP ),
387 $ KTYPE( MAXTYP )
388* ..
389* .. External Functions ..
390 DOUBLE PRECISION DLAMCH, DLARND, DSXT1
391 EXTERNAL DLAMCH, DLARND, DSXT1
392* ..
393* .. External Subroutines ..
394 EXTERNAL alasvm, dlafts, xerbla, zhbev, zhbevd,
400* ..
401* .. Intrinsic Functions ..
402 INTRINSIC abs, dble, int, log, max, min, sqrt
403* ..
404* .. Data statements ..
405 DATA ktype / 1, 2, 5*4, 5*5, 3*8, 3*9 /
406 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
407 $ 2, 3, 1, 2, 3 /
408 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
409 $ 0, 0, 4, 4, 4 /
410* ..
411* .. Executable Statements ..
412*
413* 1) Check for errors
414*
415 ntestt = 0
416 info = 0
417*
418 badnn = .false.
419 nmax = 1
420 DO 10 j = 1, nsizes
421 nmax = max( nmax, nn( j ) )
422 IF( nn( j ).LT.0 )
423 $ badnn = .true.
424 10 CONTINUE
425*
426* Check for errors
427*
428 IF( nsizes.LT.0 ) THEN
429 info = -1
430 ELSE IF( badnn ) THEN
431 info = -2
432 ELSE IF( ntypes.LT.0 ) THEN
433 info = -3
434 ELSE IF( lda.LT.nmax ) THEN
435 info = -9
436 ELSE IF( ldu.LT.nmax ) THEN
437 info = -16
438 ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
439 info = -22
440 END IF
441*
442 IF( info.NE.0 ) THEN
443 CALL xerbla( 'ZDRVST2STG', -info )
444 RETURN
445 END IF
446*
447* Quick return if nothing to do
448*
449 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
450 $ RETURN
451*
452* More Important constants
453*
454 unfl = dlamch( 'Safe minimum' )
455 ovfl = dlamch( 'Overflow' )
456 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
457 ulpinv = one / ulp
458 rtunfl = sqrt( unfl )
459 rtovfl = sqrt( ovfl )
460*
461* Loop over sizes, types
462*
463 DO 20 i = 1, 4
464 iseed2( i ) = iseed( i )
465 iseed3( i ) = iseed( i )
466 20 CONTINUE
467*
468 nerrs = 0
469 nmats = 0
470*
471 DO 1220 jsize = 1, nsizes
472 n = nn( jsize )
473 IF( n.GT.0 ) THEN
474 lgn = int( log( dble( n ) ) / log( two ) )
475 IF( 2**lgn.LT.n )
476 $ lgn = lgn + 1
477 IF( 2**lgn.LT.n )
478 $ lgn = lgn + 1
479 lwedc = max( 2*n+n*n, 2*n*n )
480 lrwedc = 1 + 4*n + 2*n*lgn + 3*n**2
481 liwedc = 3 + 5*n
482 ELSE
483 lwedc = 2
484 lrwedc = 8
485 liwedc = 8
486 END IF
487 aninv = one / dble( max( 1, n ) )
488*
489 IF( nsizes.NE.1 ) THEN
490 mtypes = min( maxtyp, ntypes )
491 ELSE
492 mtypes = min( maxtyp+1, ntypes )
493 END IF
494*
495 DO 1210 jtype = 1, mtypes
496 IF( .NOT.dotype( jtype ) )
497 $ GO TO 1210
498 nmats = nmats + 1
499 ntest = 0
500*
501 DO 30 j = 1, 4
502 ioldsd( j ) = iseed( j )
503 30 CONTINUE
504*
505* 2) Compute "A"
506*
507* Control parameters:
508*
509* KMAGN KMODE KTYPE
510* =1 O(1) clustered 1 zero
511* =2 large clustered 2 identity
512* =3 small exponential (none)
513* =4 arithmetic diagonal, (w/ eigenvalues)
514* =5 random log Hermitian, w/ eigenvalues
515* =6 random (none)
516* =7 random diagonal
517* =8 random Hermitian
518* =9 band Hermitian, w/ eigenvalues
519*
520 IF( mtypes.GT.maxtyp )
521 $ GO TO 110
522*
523 itype = ktype( jtype )
524 imode = kmode( jtype )
525*
526* Compute norm
527*
528 GO TO ( 40, 50, 60 )kmagn( jtype )
529*
530 40 CONTINUE
531 anorm = one
532 GO TO 70
533*
534 50 CONTINUE
535 anorm = ( rtovfl*ulp )*aninv
536 GO TO 70
537*
538 60 CONTINUE
539 anorm = rtunfl*n*ulpinv
540 GO TO 70
541*
542 70 CONTINUE
543*
544 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
545 iinfo = 0
546 cond = ulpinv
547*
548* Special Matrices -- Identity & Jordan block
549*
550* Zero
551*
552 IF( itype.EQ.1 ) THEN
553 iinfo = 0
554*
555 ELSE IF( itype.EQ.2 ) THEN
556*
557* Identity
558*
559 DO 80 jcol = 1, n
560 a( jcol, jcol ) = anorm
561 80 CONTINUE
562*
563 ELSE IF( itype.EQ.4 ) THEN
564*
565* Diagonal Matrix, [Eigen]values Specified
566*
567 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
568 $ anorm, 0, 0, 'N', a, lda, work, iinfo )
569*
570 ELSE IF( itype.EQ.5 ) THEN
571*
572* Hermitian, eigenvalues specified
573*
574 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
575 $ anorm, n, n, 'N', a, lda, work, iinfo )
576*
577 ELSE IF( itype.EQ.7 ) THEN
578*
579* Diagonal, random eigenvalues
580*
581 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
582 $ 'T', 'N', work( n+1 ), 1, one,
583 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
584 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
585*
586 ELSE IF( itype.EQ.8 ) THEN
587*
588* Hermitian, random eigenvalues
589*
590 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
591 $ 'T', 'N', work( n+1 ), 1, one,
592 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
593 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
594*
595 ELSE IF( itype.EQ.9 ) THEN
596*
597* Hermitian banded, eigenvalues specified
598*
599 ihbw = int( ( n-1 )*dlarnd( 1, iseed3 ) )
600 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
601 $ anorm, ihbw, ihbw, 'Z', u, ldu, work,
602 $ iinfo )
603*
604* Store as dense matrix for most routines.
605*
606 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
607 DO 100 idiag = -ihbw, ihbw
608 irow = ihbw - idiag + 1
609 j1 = max( 1, idiag+1 )
610 j2 = min( n, n+idiag )
611 DO 90 j = j1, j2
612 i = j - idiag
613 a( i, j ) = u( irow, j )
614 90 CONTINUE
615 100 CONTINUE
616 ELSE
617 iinfo = 1
618 END IF
619*
620 IF( iinfo.NE.0 ) THEN
621 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
622 $ ioldsd
623 info = abs( iinfo )
624 RETURN
625 END IF
626*
627 110 CONTINUE
628*
629 abstol = unfl + unfl
630 IF( n.LE.1 ) THEN
631 il = 1
632 iu = n
633 ELSE
634 il = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
635 iu = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
636 IF( il.GT.iu ) THEN
637 itemp = il
638 il = iu
639 iu = itemp
640 END IF
641 END IF
642*
643* Perform tests storing upper or lower triangular
644* part of matrix.
645*
646 DO 1200 iuplo = 0, 1
647 IF( iuplo.EQ.0 ) THEN
648 uplo = 'L'
649 ELSE
650 uplo = 'U'
651 END IF
652*
653* Call ZHEEVD and CHEEVX.
654*
655 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
656*
657 ntest = ntest + 1
658 CALL zheevd( 'V', uplo, n, a, ldu, d1, work, lwedc,
659 $ rwork, lrwedc, iwork, liwedc, iinfo )
660 IF( iinfo.NE.0 ) THEN
661 WRITE( nounit, fmt = 9999 )'ZHEEVD(V,' // uplo //
662 $ ')', iinfo, n, jtype, ioldsd
663 info = abs( iinfo )
664 IF( iinfo.LT.0 ) THEN
665 RETURN
666 ELSE
667 result( ntest ) = ulpinv
668 result( ntest+1 ) = ulpinv
669 result( ntest+2 ) = ulpinv
670 GO TO 130
671 END IF
672 END IF
673*
674* Do tests 1 and 2.
675*
676 CALL zhet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
677 $ ldu, tau, work, rwork, result( ntest ) )
678*
679 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
680*
681 ntest = ntest + 2
682 CALL zheevd_2stage( 'N', uplo, n, a, ldu, d3, work,
683 $ lwork, rwork, lrwedc, iwork, liwedc, iinfo )
684 IF( iinfo.NE.0 ) THEN
685 WRITE( nounit, fmt = 9999 )
686 $ 'ZHEEVD_2STAGE(N,' // uplo //
687 $ ')', iinfo, n, jtype, ioldsd
688 info = abs( iinfo )
689 IF( iinfo.LT.0 ) THEN
690 RETURN
691 ELSE
692 result( ntest ) = ulpinv
693 GO TO 130
694 END IF
695 END IF
696*
697* Do test 3.
698*
699 temp1 = zero
700 temp2 = zero
701 DO 120 j = 1, n
702 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
703 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
704 120 CONTINUE
705 result( ntest ) = temp2 / max( unfl,
706 $ ulp*max( temp1, temp2 ) )
707*
708 130 CONTINUE
709 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
710*
711 ntest = ntest + 1
712*
713 IF( n.GT.0 ) THEN
714 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
715 IF( il.NE.1 ) THEN
716 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
717 $ ten*ulp*temp3, ten*rtunfl )
718 ELSE IF( n.GT.0 ) THEN
719 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
720 $ ten*ulp*temp3, ten*rtunfl )
721 END IF
722 IF( iu.NE.n ) THEN
723 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
724 $ ten*ulp*temp3, ten*rtunfl )
725 ELSE IF( n.GT.0 ) THEN
726 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
727 $ ten*ulp*temp3, ten*rtunfl )
728 END IF
729 ELSE
730 temp3 = zero
731 vl = zero
732 vu = one
733 END IF
734*
735 CALL zheevx( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
736 $ abstol, m, wa1, z, ldu, work, lwork, rwork,
737 $ iwork, iwork( 5*n+1 ), iinfo )
738 IF( iinfo.NE.0 ) THEN
739 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,A,' // uplo //
740 $ ')', iinfo, n, jtype, ioldsd
741 info = abs( iinfo )
742 IF( iinfo.LT.0 ) THEN
743 RETURN
744 ELSE
745 result( ntest ) = ulpinv
746 result( ntest+1 ) = ulpinv
747 result( ntest+2 ) = ulpinv
748 GO TO 150
749 END IF
750 END IF
751*
752* Do tests 4 and 5.
753*
754 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
755*
756 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
757 $ ldu, tau, work, rwork, result( ntest ) )
758*
759 ntest = ntest + 2
760 CALL zheevx_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
761 $ il, iu, abstol, m2, wa2, z, ldu,
762 $ work, lwork, rwork, iwork,
763 $ iwork( 5*n+1 ), iinfo )
764 IF( iinfo.NE.0 ) THEN
765 WRITE( nounit, fmt = 9999 )
766 $ 'ZHEEVX_2STAGE(N,A,' // uplo //
767 $ ')', iinfo, n, jtype, ioldsd
768 info = abs( iinfo )
769 IF( iinfo.LT.0 ) THEN
770 RETURN
771 ELSE
772 result( ntest ) = ulpinv
773 GO TO 150
774 END IF
775 END IF
776*
777* Do test 6.
778*
779 temp1 = zero
780 temp2 = zero
781 DO 140 j = 1, n
782 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
783 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
784 140 CONTINUE
785 result( ntest ) = temp2 / max( unfl,
786 $ ulp*max( temp1, temp2 ) )
787*
788 150 CONTINUE
789 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
790*
791 ntest = ntest + 1
792*
793 CALL zheevx( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
794 $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
795 $ iwork, iwork( 5*n+1 ), iinfo )
796 IF( iinfo.NE.0 ) THEN
797 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,I,' // uplo //
798 $ ')', iinfo, n, jtype, ioldsd
799 info = abs( iinfo )
800 IF( iinfo.LT.0 ) THEN
801 RETURN
802 ELSE
803 result( ntest ) = ulpinv
804 GO TO 160
805 END IF
806 END IF
807*
808* Do tests 7 and 8.
809*
810 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
811*
812 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
813 $ v, ldu, tau, work, rwork, result( ntest ) )
814*
815 ntest = ntest + 2
816*
817 CALL zheevx_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
818 $ il, iu, abstol, m3, wa3, z, ldu,
819 $ work, lwork, rwork, iwork,
820 $ iwork( 5*n+1 ), iinfo )
821 IF( iinfo.NE.0 ) THEN
822 WRITE( nounit, fmt = 9999 )
823 $ 'ZHEEVX_2STAGE(N,I,' // uplo //
824 $ ')', iinfo, n, jtype, ioldsd
825 info = abs( iinfo )
826 IF( iinfo.LT.0 ) THEN
827 RETURN
828 ELSE
829 result( ntest ) = ulpinv
830 GO TO 160
831 END IF
832 END IF
833*
834* Do test 9.
835*
836 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
837 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
838 IF( n.GT.0 ) THEN
839 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
840 ELSE
841 temp3 = zero
842 END IF
843 result( ntest ) = ( temp1+temp2 ) /
844 $ max( unfl, temp3*ulp )
845*
846 160 CONTINUE
847 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
848*
849 ntest = ntest + 1
850*
851 CALL zheevx( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
852 $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
853 $ iwork, iwork( 5*n+1 ), iinfo )
854 IF( iinfo.NE.0 ) THEN
855 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,V,' // uplo //
856 $ ')', iinfo, n, jtype, ioldsd
857 info = abs( iinfo )
858 IF( iinfo.LT.0 ) THEN
859 RETURN
860 ELSE
861 result( ntest ) = ulpinv
862 GO TO 170
863 END IF
864 END IF
865*
866* Do tests 10 and 11.
867*
868 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
869*
870 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
871 $ v, ldu, tau, work, rwork, result( ntest ) )
872*
873 ntest = ntest + 2
874*
875 CALL zheevx_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
876 $ il, iu, abstol, m3, wa3, z, ldu,
877 $ work, lwork, rwork, iwork,
878 $ iwork( 5*n+1 ), iinfo )
879 IF( iinfo.NE.0 ) THEN
880 WRITE( nounit, fmt = 9999 )
881 $ 'ZHEEVX_2STAGE(N,V,' // uplo //
882 $ ')', iinfo, n, jtype, ioldsd
883 info = abs( iinfo )
884 IF( iinfo.LT.0 ) THEN
885 RETURN
886 ELSE
887 result( ntest ) = ulpinv
888 GO TO 170
889 END IF
890 END IF
891*
892 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
893 result( ntest ) = ulpinv
894 GO TO 170
895 END IF
896*
897* Do test 12.
898*
899 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
900 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
901 IF( n.GT.0 ) THEN
902 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
903 ELSE
904 temp3 = zero
905 END IF
906 result( ntest ) = ( temp1+temp2 ) /
907 $ max( unfl, temp3*ulp )
908*
909 170 CONTINUE
910*
911* Call ZHPEVD and CHPEVX.
912*
913 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
914*
915* Load array WORK with the upper or lower triangular
916* part of the matrix in packed form.
917*
918 IF( iuplo.EQ.1 ) THEN
919 indx = 1
920 DO 190 j = 1, n
921 DO 180 i = 1, j
922 work( indx ) = a( i, j )
923 indx = indx + 1
924 180 CONTINUE
925 190 CONTINUE
926 ELSE
927 indx = 1
928 DO 210 j = 1, n
929 DO 200 i = j, n
930 work( indx ) = a( i, j )
931 indx = indx + 1
932 200 CONTINUE
933 210 CONTINUE
934 END IF
935*
936 ntest = ntest + 1
937 indwrk = n*( n+1 ) / 2 + 1
938 CALL zhpevd( 'V', uplo, n, work, d1, z, ldu,
939 $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
940 $ liwedc, iinfo )
941 IF( iinfo.NE.0 ) THEN
942 WRITE( nounit, fmt = 9999 )'ZHPEVD(V,' // uplo //
943 $ ')', iinfo, n, jtype, ioldsd
944 info = abs( iinfo )
945 IF( iinfo.LT.0 ) THEN
946 RETURN
947 ELSE
948 result( ntest ) = ulpinv
949 result( ntest+1 ) = ulpinv
950 result( ntest+2 ) = ulpinv
951 GO TO 270
952 END IF
953 END IF
954*
955* Do tests 13 and 14.
956*
957 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
958 $ ldu, tau, work, rwork, result( ntest ) )
959*
960 IF( iuplo.EQ.1 ) THEN
961 indx = 1
962 DO 230 j = 1, n
963 DO 220 i = 1, j
964 work( indx ) = a( i, j )
965 indx = indx + 1
966 220 CONTINUE
967 230 CONTINUE
968 ELSE
969 indx = 1
970 DO 250 j = 1, n
971 DO 240 i = j, n
972 work( indx ) = a( i, j )
973 indx = indx + 1
974 240 CONTINUE
975 250 CONTINUE
976 END IF
977*
978 ntest = ntest + 2
979 indwrk = n*( n+1 ) / 2 + 1
980 CALL zhpevd( 'N', uplo, n, work, d3, z, ldu,
981 $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
982 $ liwedc, iinfo )
983 IF( iinfo.NE.0 ) THEN
984 WRITE( nounit, fmt = 9999 )'ZHPEVD(N,' // uplo //
985 $ ')', iinfo, n, jtype, ioldsd
986 info = abs( iinfo )
987 IF( iinfo.LT.0 ) THEN
988 RETURN
989 ELSE
990 result( ntest ) = ulpinv
991 GO TO 270
992 END IF
993 END IF
994*
995* Do test 15.
996*
997 temp1 = zero
998 temp2 = zero
999 DO 260 j = 1, n
1000 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1001 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1002 260 CONTINUE
1003 result( ntest ) = temp2 / max( unfl,
1004 $ ulp*max( temp1, temp2 ) )
1005*
1006* Load array WORK with the upper or lower triangular part
1007* of the matrix in packed form.
1008*
1009 270 CONTINUE
1010 IF( iuplo.EQ.1 ) THEN
1011 indx = 1
1012 DO 290 j = 1, n
1013 DO 280 i = 1, j
1014 work( indx ) = a( i, j )
1015 indx = indx + 1
1016 280 CONTINUE
1017 290 CONTINUE
1018 ELSE
1019 indx = 1
1020 DO 310 j = 1, n
1021 DO 300 i = j, n
1022 work( indx ) = a( i, j )
1023 indx = indx + 1
1024 300 CONTINUE
1025 310 CONTINUE
1026 END IF
1027*
1028 ntest = ntest + 1
1029*
1030 IF( n.GT.0 ) THEN
1031 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1032 IF( il.NE.1 ) THEN
1033 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1034 $ ten*ulp*temp3, ten*rtunfl )
1035 ELSE IF( n.GT.0 ) THEN
1036 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1037 $ ten*ulp*temp3, ten*rtunfl )
1038 END IF
1039 IF( iu.NE.n ) THEN
1040 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1041 $ ten*ulp*temp3, ten*rtunfl )
1042 ELSE IF( n.GT.0 ) THEN
1043 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1044 $ ten*ulp*temp3, ten*rtunfl )
1045 END IF
1046 ELSE
1047 temp3 = zero
1048 vl = zero
1049 vu = one
1050 END IF
1051*
1052 CALL zhpevx( 'V', 'A', uplo, n, work, vl, vu, il, iu,
1053 $ abstol, m, wa1, z, ldu, v, rwork, iwork,
1054 $ iwork( 5*n+1 ), iinfo )
1055 IF( iinfo.NE.0 ) THEN
1056 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,A,' // uplo //
1057 $ ')', iinfo, n, jtype, ioldsd
1058 info = abs( iinfo )
1059 IF( iinfo.LT.0 ) THEN
1060 RETURN
1061 ELSE
1062 result( ntest ) = ulpinv
1063 result( ntest+1 ) = ulpinv
1064 result( ntest+2 ) = ulpinv
1065 GO TO 370
1066 END IF
1067 END IF
1068*
1069* Do tests 16 and 17.
1070*
1071 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1072 $ ldu, tau, work, rwork, result( ntest ) )
1073*
1074 ntest = ntest + 2
1075*
1076 IF( iuplo.EQ.1 ) THEN
1077 indx = 1
1078 DO 330 j = 1, n
1079 DO 320 i = 1, j
1080 work( indx ) = a( i, j )
1081 indx = indx + 1
1082 320 CONTINUE
1083 330 CONTINUE
1084 ELSE
1085 indx = 1
1086 DO 350 j = 1, n
1087 DO 340 i = j, n
1088 work( indx ) = a( i, j )
1089 indx = indx + 1
1090 340 CONTINUE
1091 350 CONTINUE
1092 END IF
1093*
1094 CALL zhpevx( 'N', 'A', uplo, n, work, vl, vu, il, iu,
1095 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1096 $ iwork( 5*n+1 ), iinfo )
1097 IF( iinfo.NE.0 ) THEN
1098 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,A,' // uplo //
1099 $ ')', iinfo, n, jtype, ioldsd
1100 info = abs( iinfo )
1101 IF( iinfo.LT.0 ) THEN
1102 RETURN
1103 ELSE
1104 result( ntest ) = ulpinv
1105 GO TO 370
1106 END IF
1107 END IF
1108*
1109* Do test 18.
1110*
1111 temp1 = zero
1112 temp2 = zero
1113 DO 360 j = 1, n
1114 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1115 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1116 360 CONTINUE
1117 result( ntest ) = temp2 / max( unfl,
1118 $ ulp*max( temp1, temp2 ) )
1119*
1120 370 CONTINUE
1121 ntest = ntest + 1
1122 IF( iuplo.EQ.1 ) THEN
1123 indx = 1
1124 DO 390 j = 1, n
1125 DO 380 i = 1, j
1126 work( indx ) = a( i, j )
1127 indx = indx + 1
1128 380 CONTINUE
1129 390 CONTINUE
1130 ELSE
1131 indx = 1
1132 DO 410 j = 1, n
1133 DO 400 i = j, n
1134 work( indx ) = a( i, j )
1135 indx = indx + 1
1136 400 CONTINUE
1137 410 CONTINUE
1138 END IF
1139*
1140 CALL zhpevx( 'V', 'I', uplo, n, work, vl, vu, il, iu,
1141 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1142 $ iwork( 5*n+1 ), iinfo )
1143 IF( iinfo.NE.0 ) THEN
1144 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,I,' // uplo //
1145 $ ')', iinfo, n, jtype, ioldsd
1146 info = abs( iinfo )
1147 IF( iinfo.LT.0 ) THEN
1148 RETURN
1149 ELSE
1150 result( ntest ) = ulpinv
1151 result( ntest+1 ) = ulpinv
1152 result( ntest+2 ) = ulpinv
1153 GO TO 460
1154 END IF
1155 END IF
1156*
1157* Do tests 19 and 20.
1158*
1159 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1160 $ v, ldu, tau, work, rwork, result( ntest ) )
1161*
1162 ntest = ntest + 2
1163*
1164 IF( iuplo.EQ.1 ) THEN
1165 indx = 1
1166 DO 430 j = 1, n
1167 DO 420 i = 1, j
1168 work( indx ) = a( i, j )
1169 indx = indx + 1
1170 420 CONTINUE
1171 430 CONTINUE
1172 ELSE
1173 indx = 1
1174 DO 450 j = 1, n
1175 DO 440 i = j, n
1176 work( indx ) = a( i, j )
1177 indx = indx + 1
1178 440 CONTINUE
1179 450 CONTINUE
1180 END IF
1181*
1182 CALL zhpevx( 'N', 'I', uplo, n, work, vl, vu, il, iu,
1183 $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1184 $ iwork( 5*n+1 ), iinfo )
1185 IF( iinfo.NE.0 ) THEN
1186 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,I,' // uplo //
1187 $ ')', iinfo, n, jtype, ioldsd
1188 info = abs( iinfo )
1189 IF( iinfo.LT.0 ) THEN
1190 RETURN
1191 ELSE
1192 result( ntest ) = ulpinv
1193 GO TO 460
1194 END IF
1195 END IF
1196*
1197* Do test 21.
1198*
1199 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1200 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1201 IF( n.GT.0 ) THEN
1202 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1203 ELSE
1204 temp3 = zero
1205 END IF
1206 result( ntest ) = ( temp1+temp2 ) /
1207 $ max( unfl, temp3*ulp )
1208*
1209 460 CONTINUE
1210 ntest = ntest + 1
1211 IF( iuplo.EQ.1 ) THEN
1212 indx = 1
1213 DO 480 j = 1, n
1214 DO 470 i = 1, j
1215 work( indx ) = a( i, j )
1216 indx = indx + 1
1217 470 CONTINUE
1218 480 CONTINUE
1219 ELSE
1220 indx = 1
1221 DO 500 j = 1, n
1222 DO 490 i = j, n
1223 work( indx ) = a( i, j )
1224 indx = indx + 1
1225 490 CONTINUE
1226 500 CONTINUE
1227 END IF
1228*
1229 CALL zhpevx( 'V', 'V', uplo, n, work, vl, vu, il, iu,
1230 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1231 $ iwork( 5*n+1 ), iinfo )
1232 IF( iinfo.NE.0 ) THEN
1233 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,V,' // uplo //
1234 $ ')', iinfo, n, jtype, ioldsd
1235 info = abs( iinfo )
1236 IF( iinfo.LT.0 ) THEN
1237 RETURN
1238 ELSE
1239 result( ntest ) = ulpinv
1240 result( ntest+1 ) = ulpinv
1241 result( ntest+2 ) = ulpinv
1242 GO TO 550
1243 END IF
1244 END IF
1245*
1246* Do tests 22 and 23.
1247*
1248 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1249 $ v, ldu, tau, work, rwork, result( ntest ) )
1250*
1251 ntest = ntest + 2
1252*
1253 IF( iuplo.EQ.1 ) THEN
1254 indx = 1
1255 DO 520 j = 1, n
1256 DO 510 i = 1, j
1257 work( indx ) = a( i, j )
1258 indx = indx + 1
1259 510 CONTINUE
1260 520 CONTINUE
1261 ELSE
1262 indx = 1
1263 DO 540 j = 1, n
1264 DO 530 i = j, n
1265 work( indx ) = a( i, j )
1266 indx = indx + 1
1267 530 CONTINUE
1268 540 CONTINUE
1269 END IF
1270*
1271 CALL zhpevx( 'N', 'V', uplo, n, work, vl, vu, il, iu,
1272 $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1273 $ iwork( 5*n+1 ), iinfo )
1274 IF( iinfo.NE.0 ) THEN
1275 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,V,' // uplo //
1276 $ ')', iinfo, n, jtype, ioldsd
1277 info = abs( iinfo )
1278 IF( iinfo.LT.0 ) THEN
1279 RETURN
1280 ELSE
1281 result( ntest ) = ulpinv
1282 GO TO 550
1283 END IF
1284 END IF
1285*
1286 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1287 result( ntest ) = ulpinv
1288 GO TO 550
1289 END IF
1290*
1291* Do test 24.
1292*
1293 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1294 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1295 IF( n.GT.0 ) THEN
1296 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1297 ELSE
1298 temp3 = zero
1299 END IF
1300 result( ntest ) = ( temp1+temp2 ) /
1301 $ max( unfl, temp3*ulp )
1302*
1303 550 CONTINUE
1304*
1305* Call ZHBEVD and CHBEVX.
1306*
1307 IF( jtype.LE.7 ) THEN
1308 kd = 0
1309 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1310 kd = max( n-1, 0 )
1311 ELSE
1312 kd = ihbw
1313 END IF
1314*
1315* Load array V with the upper or lower triangular part
1316* of the matrix in band form.
1317*
1318 IF( iuplo.EQ.1 ) THEN
1319 DO 570 j = 1, n
1320 DO 560 i = max( 1, j-kd ), j
1321 v( kd+1+i-j, j ) = a( i, j )
1322 560 CONTINUE
1323 570 CONTINUE
1324 ELSE
1325 DO 590 j = 1, n
1326 DO 580 i = j, min( n, j+kd )
1327 v( 1+i-j, j ) = a( i, j )
1328 580 CONTINUE
1329 590 CONTINUE
1330 END IF
1331*
1332 ntest = ntest + 1
1333 CALL zhbevd( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1334 $ lwedc, rwork, lrwedc, iwork, liwedc, iinfo )
1335 IF( iinfo.NE.0 ) THEN
1336 WRITE( nounit, fmt = 9998 )'ZHBEVD(V,' // uplo //
1337 $ ')', iinfo, n, kd, jtype, ioldsd
1338 info = abs( iinfo )
1339 IF( iinfo.LT.0 ) THEN
1340 RETURN
1341 ELSE
1342 result( ntest ) = ulpinv
1343 result( ntest+1 ) = ulpinv
1344 result( ntest+2 ) = ulpinv
1345 GO TO 650
1346 END IF
1347 END IF
1348*
1349* Do tests 25 and 26.
1350*
1351 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1352 $ ldu, tau, work, rwork, result( ntest ) )
1353*
1354 IF( iuplo.EQ.1 ) THEN
1355 DO 610 j = 1, n
1356 DO 600 i = max( 1, j-kd ), j
1357 v( kd+1+i-j, j ) = a( i, j )
1358 600 CONTINUE
1359 610 CONTINUE
1360 ELSE
1361 DO 630 j = 1, n
1362 DO 620 i = j, min( n, j+kd )
1363 v( 1+i-j, j ) = a( i, j )
1364 620 CONTINUE
1365 630 CONTINUE
1366 END IF
1367*
1368 ntest = ntest + 2
1369 CALL zhbevd_2stage( 'N', uplo, n, kd, v, ldu, d3,
1370 $ z, ldu, work, lwork, rwork,
1371 $ lrwedc, iwork, liwedc, iinfo )
1372 IF( iinfo.NE.0 ) THEN
1373 WRITE( nounit, fmt = 9998 )
1374 $ 'ZHBEVD_2STAGE(N,' // uplo //
1375 $ ')', iinfo, n, kd, jtype, ioldsd
1376 info = abs( iinfo )
1377 IF( iinfo.LT.0 ) THEN
1378 RETURN
1379 ELSE
1380 result( ntest ) = ulpinv
1381 GO TO 650
1382 END IF
1383 END IF
1384*
1385* Do test 27.
1386*
1387 temp1 = zero
1388 temp2 = zero
1389 DO 640 j = 1, n
1390 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1391 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1392 640 CONTINUE
1393 result( ntest ) = temp2 / max( unfl,
1394 $ ulp*max( temp1, temp2 ) )
1395*
1396* Load array V with the upper or lower triangular part
1397* of the matrix in band form.
1398*
1399 650 CONTINUE
1400 IF( iuplo.EQ.1 ) THEN
1401 DO 670 j = 1, n
1402 DO 660 i = max( 1, j-kd ), j
1403 v( kd+1+i-j, j ) = a( i, j )
1404 660 CONTINUE
1405 670 CONTINUE
1406 ELSE
1407 DO 690 j = 1, n
1408 DO 680 i = j, min( n, j+kd )
1409 v( 1+i-j, j ) = a( i, j )
1410 680 CONTINUE
1411 690 CONTINUE
1412 END IF
1413*
1414 ntest = ntest + 1
1415 CALL zhbevx( 'V', 'A', uplo, n, kd, v, ldu, u, ldu, vl,
1416 $ vu, il, iu, abstol, m, wa1, z, ldu, work,
1417 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1418 IF( iinfo.NE.0 ) THEN
1419 WRITE( nounit, fmt = 9999 )'ZHBEVX(V,A,' // uplo //
1420 $ ')', iinfo, n, kd, jtype, ioldsd
1421 info = abs( iinfo )
1422 IF( iinfo.LT.0 ) THEN
1423 RETURN
1424 ELSE
1425 result( ntest ) = ulpinv
1426 result( ntest+1 ) = ulpinv
1427 result( ntest+2 ) = ulpinv
1428 GO TO 750
1429 END IF
1430 END IF
1431*
1432* Do tests 28 and 29.
1433*
1434 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1435 $ ldu, tau, work, rwork, result( ntest ) )
1436*
1437 ntest = ntest + 2
1438*
1439 IF( iuplo.EQ.1 ) THEN
1440 DO 710 j = 1, n
1441 DO 700 i = max( 1, j-kd ), j
1442 v( kd+1+i-j, j ) = a( i, j )
1443 700 CONTINUE
1444 710 CONTINUE
1445 ELSE
1446 DO 730 j = 1, n
1447 DO 720 i = j, min( n, j+kd )
1448 v( 1+i-j, j ) = a( i, j )
1449 720 CONTINUE
1450 730 CONTINUE
1451 END IF
1452*
1453 CALL zhbevx_2stage( 'N', 'A', uplo, n, kd, v, ldu,
1454 $ u, ldu, vl, vu, il, iu, abstol,
1455 $ m2, wa2, z, ldu, work, lwork,
1456 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1457 IF( iinfo.NE.0 ) THEN
1458 WRITE( nounit, fmt = 9998 )
1459 $ 'ZHBEVX_2STAGE(N,A,' // uplo //
1460 $ ')', iinfo, n, kd, jtype, ioldsd
1461 info = abs( iinfo )
1462 IF( iinfo.LT.0 ) THEN
1463 RETURN
1464 ELSE
1465 result( ntest ) = ulpinv
1466 GO TO 750
1467 END IF
1468 END IF
1469*
1470* Do test 30.
1471*
1472 temp1 = zero
1473 temp2 = zero
1474 DO 740 j = 1, n
1475 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1476 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1477 740 CONTINUE
1478 result( ntest ) = temp2 / max( unfl,
1479 $ ulp*max( temp1, temp2 ) )
1480*
1481* Load array V with the upper or lower triangular part
1482* of the matrix in band form.
1483*
1484 750 CONTINUE
1485 ntest = ntest + 1
1486 IF( iuplo.EQ.1 ) THEN
1487 DO 770 j = 1, n
1488 DO 760 i = max( 1, j-kd ), j
1489 v( kd+1+i-j, j ) = a( i, j )
1490 760 CONTINUE
1491 770 CONTINUE
1492 ELSE
1493 DO 790 j = 1, n
1494 DO 780 i = j, min( n, j+kd )
1495 v( 1+i-j, j ) = a( i, j )
1496 780 CONTINUE
1497 790 CONTINUE
1498 END IF
1499*
1500 CALL zhbevx( 'V', 'I', uplo, n, kd, v, ldu, u, ldu, vl,
1501 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1502 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1503 IF( iinfo.NE.0 ) THEN
1504 WRITE( nounit, fmt = 9998 )'ZHBEVX(V,I,' // uplo //
1505 $ ')', iinfo, n, kd, jtype, ioldsd
1506 info = abs( iinfo )
1507 IF( iinfo.LT.0 ) THEN
1508 RETURN
1509 ELSE
1510 result( ntest ) = ulpinv
1511 result( ntest+1 ) = ulpinv
1512 result( ntest+2 ) = ulpinv
1513 GO TO 840
1514 END IF
1515 END IF
1516*
1517* Do tests 31 and 32.
1518*
1519 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1520 $ v, ldu, tau, work, rwork, result( ntest ) )
1521*
1522 ntest = ntest + 2
1523*
1524 IF( iuplo.EQ.1 ) THEN
1525 DO 810 j = 1, n
1526 DO 800 i = max( 1, j-kd ), j
1527 v( kd+1+i-j, j ) = a( i, j )
1528 800 CONTINUE
1529 810 CONTINUE
1530 ELSE
1531 DO 830 j = 1, n
1532 DO 820 i = j, min( n, j+kd )
1533 v( 1+i-j, j ) = a( i, j )
1534 820 CONTINUE
1535 830 CONTINUE
1536 END IF
1537 CALL zhbevx_2stage( 'N', 'I', uplo, n, kd, v, ldu,
1538 $ u, ldu, vl, vu, il, iu, abstol,
1539 $ m3, wa3, z, ldu, work, lwork,
1540 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1541 IF( iinfo.NE.0 ) THEN
1542 WRITE( nounit, fmt = 9998 )
1543 $ 'ZHBEVX_2STAGE(N,I,' // uplo //
1544 $ ')', iinfo, n, kd, jtype, ioldsd
1545 info = abs( iinfo )
1546 IF( iinfo.LT.0 ) THEN
1547 RETURN
1548 ELSE
1549 result( ntest ) = ulpinv
1550 GO TO 840
1551 END IF
1552 END IF
1553*
1554* Do test 33.
1555*
1556 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1557 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1558 IF( n.GT.0 ) THEN
1559 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1560 ELSE
1561 temp3 = zero
1562 END IF
1563 result( ntest ) = ( temp1+temp2 ) /
1564 $ max( unfl, temp3*ulp )
1565*
1566* Load array V with the upper or lower triangular part
1567* of the matrix in band form.
1568*
1569 840 CONTINUE
1570 ntest = ntest + 1
1571 IF( iuplo.EQ.1 ) THEN
1572 DO 860 j = 1, n
1573 DO 850 i = max( 1, j-kd ), j
1574 v( kd+1+i-j, j ) = a( i, j )
1575 850 CONTINUE
1576 860 CONTINUE
1577 ELSE
1578 DO 880 j = 1, n
1579 DO 870 i = j, min( n, j+kd )
1580 v( 1+i-j, j ) = a( i, j )
1581 870 CONTINUE
1582 880 CONTINUE
1583 END IF
1584 CALL zhbevx( 'V', 'V', uplo, n, kd, v, ldu, u, ldu, vl,
1585 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1586 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1587 IF( iinfo.NE.0 ) THEN
1588 WRITE( nounit, fmt = 9998 )'ZHBEVX(V,V,' // uplo //
1589 $ ')', iinfo, n, kd, jtype, ioldsd
1590 info = abs( iinfo )
1591 IF( iinfo.LT.0 ) THEN
1592 RETURN
1593 ELSE
1594 result( ntest ) = ulpinv
1595 result( ntest+1 ) = ulpinv
1596 result( ntest+2 ) = ulpinv
1597 GO TO 930
1598 END IF
1599 END IF
1600*
1601* Do tests 34 and 35.
1602*
1603 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1604 $ v, ldu, tau, work, rwork, result( ntest ) )
1605*
1606 ntest = ntest + 2
1607*
1608 IF( iuplo.EQ.1 ) THEN
1609 DO 900 j = 1, n
1610 DO 890 i = max( 1, j-kd ), j
1611 v( kd+1+i-j, j ) = a( i, j )
1612 890 CONTINUE
1613 900 CONTINUE
1614 ELSE
1615 DO 920 j = 1, n
1616 DO 910 i = j, min( n, j+kd )
1617 v( 1+i-j, j ) = a( i, j )
1618 910 CONTINUE
1619 920 CONTINUE
1620 END IF
1621 CALL zhbevx_2stage( 'N', 'V', uplo, n, kd, v, ldu,
1622 $ u, ldu, vl, vu, il, iu, abstol,
1623 $ m3, wa3, z, ldu, work, lwork,
1624 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1625 IF( iinfo.NE.0 ) THEN
1626 WRITE( nounit, fmt = 9998 )
1627 $ 'ZHBEVX_2STAGE(N,V,' // uplo //
1628 $ ')', iinfo, n, kd, jtype, ioldsd
1629 info = abs( iinfo )
1630 IF( iinfo.LT.0 ) THEN
1631 RETURN
1632 ELSE
1633 result( ntest ) = ulpinv
1634 GO TO 930
1635 END IF
1636 END IF
1637*
1638 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1639 result( ntest ) = ulpinv
1640 GO TO 930
1641 END IF
1642*
1643* Do test 36.
1644*
1645 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1646 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1647 IF( n.GT.0 ) THEN
1648 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1649 ELSE
1650 temp3 = zero
1651 END IF
1652 result( ntest ) = ( temp1+temp2 ) /
1653 $ max( unfl, temp3*ulp )
1654*
1655 930 CONTINUE
1656*
1657* Call ZHEEV
1658*
1659 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
1660*
1661 ntest = ntest + 1
1662 CALL zheev( 'V', uplo, n, a, ldu, d1, work, lwork, rwork,
1663 $ iinfo )
1664 IF( iinfo.NE.0 ) THEN
1665 WRITE( nounit, fmt = 9999 )'ZHEEV(V,' // uplo // ')',
1666 $ iinfo, n, jtype, ioldsd
1667 info = abs( iinfo )
1668 IF( iinfo.LT.0 ) THEN
1669 RETURN
1670 ELSE
1671 result( ntest ) = ulpinv
1672 result( ntest+1 ) = ulpinv
1673 result( ntest+2 ) = ulpinv
1674 GO TO 950
1675 END IF
1676 END IF
1677*
1678* Do tests 37 and 38
1679*
1680 CALL zhet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
1681 $ ldu, tau, work, rwork, result( ntest ) )
1682*
1683 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1684*
1685 ntest = ntest + 2
1686 CALL zheev_2stage( 'N', uplo, n, a, ldu, d3,
1687 $ work, lwork, rwork, iinfo )
1688 IF( iinfo.NE.0 ) THEN
1689 WRITE( nounit, fmt = 9999 )
1690 $ 'ZHEEV_2STAGE(N,' // uplo // ')',
1691 $ iinfo, n, jtype, ioldsd
1692 info = abs( iinfo )
1693 IF( iinfo.LT.0 ) THEN
1694 RETURN
1695 ELSE
1696 result( ntest ) = ulpinv
1697 GO TO 950
1698 END IF
1699 END IF
1700*
1701* Do test 39
1702*
1703 temp1 = zero
1704 temp2 = zero
1705 DO 940 j = 1, n
1706 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1707 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1708 940 CONTINUE
1709 result( ntest ) = temp2 / max( unfl,
1710 $ ulp*max( temp1, temp2 ) )
1711*
1712 950 CONTINUE
1713*
1714 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1715*
1716* Call ZHPEV
1717*
1718* Load array WORK with the upper or lower triangular
1719* part of the matrix in packed form.
1720*
1721 IF( iuplo.EQ.1 ) THEN
1722 indx = 1
1723 DO 970 j = 1, n
1724 DO 960 i = 1, j
1725 work( indx ) = a( i, j )
1726 indx = indx + 1
1727 960 CONTINUE
1728 970 CONTINUE
1729 ELSE
1730 indx = 1
1731 DO 990 j = 1, n
1732 DO 980 i = j, n
1733 work( indx ) = a( i, j )
1734 indx = indx + 1
1735 980 CONTINUE
1736 990 CONTINUE
1737 END IF
1738*
1739 ntest = ntest + 1
1740 indwrk = n*( n+1 ) / 2 + 1
1741 CALL zhpev( 'V', uplo, n, work, d1, z, ldu,
1742 $ work( indwrk ), rwork, iinfo )
1743 IF( iinfo.NE.0 ) THEN
1744 WRITE( nounit, fmt = 9999 )'ZHPEV(V,' // uplo // ')',
1745 $ iinfo, n, jtype, ioldsd
1746 info = abs( iinfo )
1747 IF( iinfo.LT.0 ) THEN
1748 RETURN
1749 ELSE
1750 result( ntest ) = ulpinv
1751 result( ntest+1 ) = ulpinv
1752 result( ntest+2 ) = ulpinv
1753 GO TO 1050
1754 END IF
1755 END IF
1756*
1757* Do tests 40 and 41.
1758*
1759 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1760 $ ldu, tau, work, rwork, result( ntest ) )
1761*
1762 IF( iuplo.EQ.1 ) THEN
1763 indx = 1
1764 DO 1010 j = 1, n
1765 DO 1000 i = 1, j
1766 work( indx ) = a( i, j )
1767 indx = indx + 1
1768 1000 CONTINUE
1769 1010 CONTINUE
1770 ELSE
1771 indx = 1
1772 DO 1030 j = 1, n
1773 DO 1020 i = j, n
1774 work( indx ) = a( i, j )
1775 indx = indx + 1
1776 1020 CONTINUE
1777 1030 CONTINUE
1778 END IF
1779*
1780 ntest = ntest + 2
1781 indwrk = n*( n+1 ) / 2 + 1
1782 CALL zhpev( 'N', uplo, n, work, d3, z, ldu,
1783 $ work( indwrk ), rwork, iinfo )
1784 IF( iinfo.NE.0 ) THEN
1785 WRITE( nounit, fmt = 9999 )'ZHPEV(N,' // uplo // ')',
1786 $ iinfo, n, jtype, ioldsd
1787 info = abs( iinfo )
1788 IF( iinfo.LT.0 ) THEN
1789 RETURN
1790 ELSE
1791 result( ntest ) = ulpinv
1792 GO TO 1050
1793 END IF
1794 END IF
1795*
1796* Do test 42
1797*
1798 temp1 = zero
1799 temp2 = zero
1800 DO 1040 j = 1, n
1801 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1802 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1803 1040 CONTINUE
1804 result( ntest ) = temp2 / max( unfl,
1805 $ ulp*max( temp1, temp2 ) )
1806*
1807 1050 CONTINUE
1808*
1809* Call ZHBEV
1810*
1811 IF( jtype.LE.7 ) THEN
1812 kd = 0
1813 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1814 kd = max( n-1, 0 )
1815 ELSE
1816 kd = ihbw
1817 END IF
1818*
1819* Load array V with the upper or lower triangular part
1820* of the matrix in band form.
1821*
1822 IF( iuplo.EQ.1 ) THEN
1823 DO 1070 j = 1, n
1824 DO 1060 i = max( 1, j-kd ), j
1825 v( kd+1+i-j, j ) = a( i, j )
1826 1060 CONTINUE
1827 1070 CONTINUE
1828 ELSE
1829 DO 1090 j = 1, n
1830 DO 1080 i = j, min( n, j+kd )
1831 v( 1+i-j, j ) = a( i, j )
1832 1080 CONTINUE
1833 1090 CONTINUE
1834 END IF
1835*
1836 ntest = ntest + 1
1837 CALL zhbev( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1838 $ rwork, iinfo )
1839 IF( iinfo.NE.0 ) THEN
1840 WRITE( nounit, fmt = 9998 )'ZHBEV(V,' // uplo // ')',
1841 $ iinfo, n, kd, jtype, ioldsd
1842 info = abs( iinfo )
1843 IF( iinfo.LT.0 ) THEN
1844 RETURN
1845 ELSE
1846 result( ntest ) = ulpinv
1847 result( ntest+1 ) = ulpinv
1848 result( ntest+2 ) = ulpinv
1849 GO TO 1140
1850 END IF
1851 END IF
1852*
1853* Do tests 43 and 44.
1854*
1855 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1856 $ ldu, tau, work, rwork, result( ntest ) )
1857*
1858 IF( iuplo.EQ.1 ) THEN
1859 DO 1110 j = 1, n
1860 DO 1100 i = max( 1, j-kd ), j
1861 v( kd+1+i-j, j ) = a( i, j )
1862 1100 CONTINUE
1863 1110 CONTINUE
1864 ELSE
1865 DO 1130 j = 1, n
1866 DO 1120 i = j, min( n, j+kd )
1867 v( 1+i-j, j ) = a( i, j )
1868 1120 CONTINUE
1869 1130 CONTINUE
1870 END IF
1871*
1872 ntest = ntest + 2
1873 CALL zhbev_2stage( 'N', uplo, n, kd, v, ldu, d3, z, ldu,
1874 $ work, lwork, rwork, iinfo )
1875 IF( iinfo.NE.0 ) THEN
1876 WRITE( nounit, fmt = 9998 )
1877 $ 'ZHBEV_2STAGE(N,' // uplo // ')',
1878 $ iinfo, n, kd, jtype, ioldsd
1879 info = abs( iinfo )
1880 IF( iinfo.LT.0 ) THEN
1881 RETURN
1882 ELSE
1883 result( ntest ) = ulpinv
1884 GO TO 1140
1885 END IF
1886 END IF
1887*
1888 1140 CONTINUE
1889*
1890* Do test 45.
1891*
1892 temp1 = zero
1893 temp2 = zero
1894 DO 1150 j = 1, n
1895 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1896 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1897 1150 CONTINUE
1898 result( ntest ) = temp2 / max( unfl,
1899 $ ulp*max( temp1, temp2 ) )
1900*
1901 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
1902 ntest = ntest + 1
1903 CALL zheevr( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
1904 $ abstol, m, wa1, z, ldu, iwork, work, lwork,
1905 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1906 $ iinfo )
1907 IF( iinfo.NE.0 ) THEN
1908 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,A,' // uplo //
1909 $ ')', iinfo, n, jtype, ioldsd
1910 info = abs( iinfo )
1911 IF( iinfo.LT.0 ) THEN
1912 RETURN
1913 ELSE
1914 result( ntest ) = ulpinv
1915 result( ntest+1 ) = ulpinv
1916 result( ntest+2 ) = ulpinv
1917 GO TO 1170
1918 END IF
1919 END IF
1920*
1921* Do tests 45 and 46 (or ... )
1922*
1923 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1924*
1925 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1926 $ ldu, tau, work, rwork, result( ntest ) )
1927*
1928 ntest = ntest + 2
1929 CALL zheevr_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
1930 $ il, iu, abstol, m2, wa2, z, ldu,
1931 $ iwork, work, lwork, rwork, lrwork,
1932 $ iwork( 2*n+1 ), liwork-2*n, iinfo )
1933 IF( iinfo.NE.0 ) THEN
1934 WRITE( nounit, fmt = 9999 )
1935 $ 'ZHEEVR_2STAGE(N,A,' // uplo //
1936 $ ')', iinfo, n, jtype, ioldsd
1937 info = abs( iinfo )
1938 IF( iinfo.LT.0 ) THEN
1939 RETURN
1940 ELSE
1941 result( ntest ) = ulpinv
1942 GO TO 1170
1943 END IF
1944 END IF
1945*
1946* Do test 47 (or ... )
1947*
1948 temp1 = zero
1949 temp2 = zero
1950 DO 1160 j = 1, n
1951 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1952 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1953 1160 CONTINUE
1954 result( ntest ) = temp2 / max( unfl,
1955 $ ulp*max( temp1, temp2 ) )
1956*
1957 1170 CONTINUE
1958*
1959 ntest = ntest + 1
1960 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1961 CALL zheevr( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
1962 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
1963 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1964 $ iinfo )
1965 IF( iinfo.NE.0 ) THEN
1966 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,I,' // uplo //
1967 $ ')', iinfo, n, jtype, ioldsd
1968 info = abs( iinfo )
1969 IF( iinfo.LT.0 ) THEN
1970 RETURN
1971 ELSE
1972 result( ntest ) = ulpinv
1973 result( ntest+1 ) = ulpinv
1974 result( ntest+2 ) = ulpinv
1975 GO TO 1180
1976 END IF
1977 END IF
1978*
1979* Do tests 48 and 49 (or +??)
1980*
1981 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1982*
1983 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1984 $ v, ldu, tau, work, rwork, result( ntest ) )
1985*
1986 ntest = ntest + 2
1987 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1988 CALL zheevr_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
1989 $ il, iu, abstol, m3, wa3, z, ldu,
1990 $ iwork, work, lwork, rwork, lrwork,
1991 $ iwork( 2*n+1 ), liwork-2*n, iinfo )
1992 IF( iinfo.NE.0 ) THEN
1993 WRITE( nounit, fmt = 9999 )
1994 $ 'ZHEEVR_2STAGE(N,I,' // uplo //
1995 $ ')', iinfo, n, jtype, ioldsd
1996 info = abs( iinfo )
1997 IF( iinfo.LT.0 ) THEN
1998 RETURN
1999 ELSE
2000 result( ntest ) = ulpinv
2001 GO TO 1180
2002 END IF
2003 END IF
2004*
2005* Do test 50 (or +??)
2006*
2007 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2008 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2009 result( ntest ) = ( temp1+temp2 ) /
2010 $ max( unfl, ulp*temp3 )
2011 1180 CONTINUE
2012*
2013 ntest = ntest + 1
2014 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2015 CALL zheevr( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
2016 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2017 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
2018 $ iinfo )
2019 IF( iinfo.NE.0 ) THEN
2020 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,V,' // uplo //
2021 $ ')', iinfo, n, jtype, ioldsd
2022 info = abs( iinfo )
2023 IF( iinfo.LT.0 ) THEN
2024 RETURN
2025 ELSE
2026 result( ntest ) = ulpinv
2027 result( ntest+1 ) = ulpinv
2028 result( ntest+2 ) = ulpinv
2029 GO TO 1190
2030 END IF
2031 END IF
2032*
2033* Do tests 51 and 52 (or +??)
2034*
2035 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2036*
2037 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2038 $ v, ldu, tau, work, rwork, result( ntest ) )
2039*
2040 ntest = ntest + 2
2041 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2042 CALL zheevr_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
2043 $ il, iu, abstol, m3, wa3, z, ldu,
2044 $ iwork, work, lwork, rwork, lrwork,
2045 $ iwork( 2*n+1 ), liwork-2*n, iinfo )
2046 IF( iinfo.NE.0 ) THEN
2047 WRITE( nounit, fmt = 9999 )
2048 $ 'ZHEEVR_2STAGE(N,V,' // uplo //
2049 $ ')', iinfo, n, jtype, ioldsd
2050 info = abs( iinfo )
2051 IF( iinfo.LT.0 ) THEN
2052 RETURN
2053 ELSE
2054 result( ntest ) = ulpinv
2055 GO TO 1190
2056 END IF
2057 END IF
2058*
2059 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2060 result( ntest ) = ulpinv
2061 GO TO 1190
2062 END IF
2063*
2064* Do test 52 (or +??)
2065*
2066 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2067 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2068 IF( n.GT.0 ) THEN
2069 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2070 ELSE
2071 temp3 = zero
2072 END IF
2073 result( ntest ) = ( temp1+temp2 ) /
2074 $ max( unfl, temp3*ulp )
2075*
2076 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2077*
2078*
2079*
2080*
2081* Load array V with the upper or lower triangular part
2082* of the matrix in band form.
2083*
2084 1190 CONTINUE
2085*
2086 1200 CONTINUE
2087*
2088* End of Loop -- Check for RESULT(j) > THRESH
2089*
2090 ntestt = ntestt + ntest
2091 CALL dlafts( 'ZST', n, n, jtype, ntest, result, ioldsd,
2092 $ thresh, nounit, nerrs )
2093*
2094 1210 CONTINUE
2095 1220 CONTINUE
2096*
2097* Summary
2098*
2099 CALL alasvm( 'ZST', nounit, nerrs, ntestt, 0 )
2100*
2101 9999 FORMAT( ' ZDRVST2STG: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2102 $ ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
2103 9998 FORMAT( ' ZDRVST2STG: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2104 $ ', KD=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
2105 $ ')' )
2106*
2107 RETURN
2108*
2109* End of ZDRVST2STG
2110*
2111 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
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 zhbev_2stage(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, rwork, info)
ZHBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER m...
subroutine zhbev(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, rwork, info)
ZHBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
Definition zhbev.f:152
subroutine zhbevd_2stage(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine zhbevd(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHBEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition zhbevd.f:209
subroutine zhbevx_2stage(jobz, range, uplo, n, kd, ab, ldab, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, ifail, info)
ZHBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine zhbevx(jobz, range, uplo, n, kd, ab, ldab, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
ZHBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition zhbevx.f:267
subroutine zheev_2stage(jobz, uplo, n, a, lda, w, work, lwork, rwork, info)
ZHEEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matr...
subroutine zheev(jobz, uplo, n, a, lda, w, work, lwork, rwork, info)
ZHEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheev.f:140
subroutine zheevd_2stage(jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHEEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat...
subroutine zheevd(jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHEEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheevd.f:199
subroutine zheevr_2stage(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHEEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat...
subroutine zheevr(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheevr.f:357
subroutine zheevx_2stage(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, ifail, info)
ZHEEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat...
subroutine zheevx(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, ifail, info)
ZHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheevx.f:259
subroutine zhetrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
ZHETRD_2STAGE
subroutine zhpev(jobz, uplo, n, ap, w, z, ldz, work, rwork, info)
ZHPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
Definition zhpev.f:138
subroutine zhpevd(jobz, uplo, n, ap, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition zhpevd.f:194
subroutine zhpevx(jobz, range, uplo, n, ap, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
ZHPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition zhpevx.f:240
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 zdrvst2stg(nsizes, nn, ntypes, dotype, iseed, thresh, nounit, a, lda, d1, d2, d3, wa1, wa2, wa3, u, ldu, v, tau, z, work, lwork, rwork, lrwork, iwork, liwork, result, info)
ZDRVST2STG
Definition zdrvst2stg.f:338
subroutine zhet21(itype, uplo, n, kband, a, lda, d, e, u, ldu, v, ldv, tau, work, rwork, result)
ZHET21
Definition zhet21.f:214
subroutine zhet22(itype, uplo, n, m, kband, a, lda, d, e, u, ldu, v, ldv, tau, work, rwork, result)
ZHET22
Definition zhet22.f:161
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