LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvst.f
Go to the documentation of this file.
1*> \brief \b ZDRVST
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 ZDRVST( 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*> ZDRVST 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 ZDRVST 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*> ZDRVST 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, ZDRVST
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 ZDRVST 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 zdrvst( 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, zhbevx,
397 $ zlatms
398* ..
399* .. Intrinsic Functions ..
400 INTRINSIC abs, dble, int, log, max, min, sqrt
401* ..
402* .. Data statements ..
403 DATA ktype / 1, 2, 5*4, 5*5, 3*8, 3*9 /
404 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
405 $ 2, 3, 1, 2, 3 /
406 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
407 $ 0, 0, 4, 4, 4 /
408* ..
409* .. Executable Statements ..
410*
411* 1) Check for errors
412*
413 ntestt = 0
414 info = 0
415*
416 badnn = .false.
417 nmax = 1
418 DO 10 j = 1, nsizes
419 nmax = max( nmax, nn( j ) )
420 IF( nn( j ).LT.0 )
421 $ badnn = .true.
422 10 CONTINUE
423*
424* Check for errors
425*
426 IF( nsizes.LT.0 ) THEN
427 info = -1
428 ELSE IF( badnn ) THEN
429 info = -2
430 ELSE IF( ntypes.LT.0 ) THEN
431 info = -3
432 ELSE IF( lda.LT.nmax ) THEN
433 info = -9
434 ELSE IF( ldu.LT.nmax ) THEN
435 info = -16
436 ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
437 info = -22
438 END IF
439*
440 IF( info.NE.0 ) THEN
441 CALL xerbla( 'ZDRVST', -info )
442 RETURN
443 END IF
444*
445* Quick return if nothing to do
446*
447 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
448 $ RETURN
449*
450* More Important constants
451*
452 unfl = dlamch( 'Safe minimum' )
453 ovfl = dlamch( 'Overflow' )
454 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
455 ulpinv = one / ulp
456 rtunfl = sqrt( unfl )
457 rtovfl = sqrt( ovfl )
458*
459* Loop over sizes, types
460*
461 DO 20 i = 1, 4
462 iseed2( i ) = iseed( i )
463 iseed3( i ) = iseed( i )
464 20 CONTINUE
465*
466 nerrs = 0
467 nmats = 0
468*
469 DO 1220 jsize = 1, nsizes
470 n = nn( jsize )
471 IF( n.GT.0 ) THEN
472 lgn = int( log( dble( n ) ) / log( two ) )
473 IF( 2**lgn.LT.n )
474 $ lgn = lgn + 1
475 IF( 2**lgn.LT.n )
476 $ lgn = lgn + 1
477 lwedc = max( 2*n+n*n, 2*n*n )
478 lrwedc = 1 + 4*n + 2*n*lgn + 3*n**2
479 liwedc = 3 + 5*n
480 ELSE
481 lwedc = 2
482 lrwedc = 8
483 liwedc = 8
484 END IF
485 aninv = one / dble( max( 1, n ) )
486*
487 IF( nsizes.NE.1 ) THEN
488 mtypes = min( maxtyp, ntypes )
489 ELSE
490 mtypes = min( maxtyp+1, ntypes )
491 END IF
492*
493 DO 1210 jtype = 1, mtypes
494 IF( .NOT.dotype( jtype ) )
495 $ GO TO 1210
496 nmats = nmats + 1
497 ntest = 0
498*
499 DO 30 j = 1, 4
500 ioldsd( j ) = iseed( j )
501 30 CONTINUE
502*
503* 2) Compute "A"
504*
505* Control parameters:
506*
507* KMAGN KMODE KTYPE
508* =1 O(1) clustered 1 zero
509* =2 large clustered 2 identity
510* =3 small exponential (none)
511* =4 arithmetic diagonal, (w/ eigenvalues)
512* =5 random log Hermitian, w/ eigenvalues
513* =6 random (none)
514* =7 random diagonal
515* =8 random Hermitian
516* =9 band Hermitian, w/ eigenvalues
517*
518 IF( mtypes.GT.maxtyp )
519 $ GO TO 110
520*
521 itype = ktype( jtype )
522 imode = kmode( jtype )
523*
524* Compute norm
525*
526 GO TO ( 40, 50, 60 )kmagn( jtype )
527*
528 40 CONTINUE
529 anorm = one
530 GO TO 70
531*
532 50 CONTINUE
533 anorm = ( rtovfl*ulp )*aninv
534 GO TO 70
535*
536 60 CONTINUE
537 anorm = rtunfl*n*ulpinv
538 GO TO 70
539*
540 70 CONTINUE
541*
542 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
543 iinfo = 0
544 cond = ulpinv
545*
546* Special Matrices -- Identity & Jordan block
547*
548* Zero
549*
550 IF( itype.EQ.1 ) THEN
551 iinfo = 0
552*
553 ELSE IF( itype.EQ.2 ) THEN
554*
555* Identity
556*
557 DO 80 jcol = 1, n
558 a( jcol, jcol ) = anorm
559 80 CONTINUE
560*
561 ELSE IF( itype.EQ.4 ) THEN
562*
563* Diagonal Matrix, [Eigen]values Specified
564*
565 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
566 $ anorm, 0, 0, 'N', a, lda, work, iinfo )
567*
568 ELSE IF( itype.EQ.5 ) THEN
569*
570* Hermitian, eigenvalues specified
571*
572 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
573 $ anorm, n, n, 'N', a, lda, work, iinfo )
574*
575 ELSE IF( itype.EQ.7 ) THEN
576*
577* Diagonal, random eigenvalues
578*
579 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
580 $ 'T', 'N', work( n+1 ), 1, one,
581 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
582 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
583*
584 ELSE IF( itype.EQ.8 ) THEN
585*
586* Hermitian, random eigenvalues
587*
588 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
589 $ 'T', 'N', work( n+1 ), 1, one,
590 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
591 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
592*
593 ELSE IF( itype.EQ.9 ) THEN
594*
595* Hermitian banded, eigenvalues specified
596*
597 ihbw = int( ( n-1 )*dlarnd( 1, iseed3 ) )
598 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
599 $ anorm, ihbw, ihbw, 'Z', u, ldu, work,
600 $ iinfo )
601*
602* Store as dense matrix for most routines.
603*
604 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
605 DO 100 idiag = -ihbw, ihbw
606 irow = ihbw - idiag + 1
607 j1 = max( 1, idiag+1 )
608 j2 = min( n, n+idiag )
609 DO 90 j = j1, j2
610 i = j - idiag
611 a( i, j ) = u( irow, j )
612 90 CONTINUE
613 100 CONTINUE
614 ELSE
615 iinfo = 1
616 END IF
617*
618 IF( iinfo.NE.0 ) THEN
619 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
620 $ ioldsd
621 info = abs( iinfo )
622 RETURN
623 END IF
624*
625 110 CONTINUE
626*
627 abstol = unfl + unfl
628 IF( n.LE.1 ) THEN
629 il = 1
630 iu = n
631 ELSE
632 il = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
633 iu = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
634 IF( il.GT.iu ) THEN
635 itemp = il
636 il = iu
637 iu = itemp
638 END IF
639 END IF
640*
641* Perform tests storing upper or lower triangular
642* part of matrix.
643*
644 DO 1200 iuplo = 0, 1
645 IF( iuplo.EQ.0 ) THEN
646 uplo = 'L'
647 ELSE
648 uplo = 'U'
649 END IF
650*
651* Call ZHEEVD and CHEEVX.
652*
653 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
654*
655 ntest = ntest + 1
656 CALL zheevd( 'V', uplo, n, a, ldu, d1, work, lwedc,
657 $ rwork, lrwedc, iwork, liwedc, iinfo )
658 IF( iinfo.NE.0 ) THEN
659 WRITE( nounit, fmt = 9999 )'ZHEEVD(V,' // uplo //
660 $ ')', iinfo, n, jtype, ioldsd
661 info = abs( iinfo )
662 IF( iinfo.LT.0 ) THEN
663 RETURN
664 ELSE
665 result( ntest ) = ulpinv
666 result( ntest+1 ) = ulpinv
667 result( ntest+2 ) = ulpinv
668 GO TO 130
669 END IF
670 END IF
671*
672* Do tests 1 and 2.
673*
674 CALL zhet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
675 $ ldu, tau, work, rwork, result( ntest ) )
676*
677 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
678*
679 ntest = ntest + 2
680 CALL zheevd( 'N', uplo, n, a, ldu, d3, work, lwedc,
681 $ rwork, lrwedc, iwork, liwedc, iinfo )
682 IF( iinfo.NE.0 ) THEN
683 WRITE( nounit, fmt = 9999 )'ZHEEVD(N,' // uplo //
684 $ ')', iinfo, n, jtype, ioldsd
685 info = abs( iinfo )
686 IF( iinfo.LT.0 ) THEN
687 RETURN
688 ELSE
689 result( ntest ) = ulpinv
690 GO TO 130
691 END IF
692 END IF
693*
694* Do test 3.
695*
696 temp1 = zero
697 temp2 = zero
698 DO 120 j = 1, n
699 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
700 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
701 120 CONTINUE
702 result( ntest ) = temp2 / max( unfl,
703 $ ulp*max( temp1, temp2 ) )
704*
705 130 CONTINUE
706 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
707*
708 ntest = ntest + 1
709*
710 IF( n.GT.0 ) THEN
711 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
712 IF( il.NE.1 ) THEN
713 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
714 $ ten*ulp*temp3, ten*rtunfl )
715 ELSE IF( n.GT.0 ) THEN
716 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
717 $ ten*ulp*temp3, ten*rtunfl )
718 END IF
719 IF( iu.NE.n ) THEN
720 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
721 $ ten*ulp*temp3, ten*rtunfl )
722 ELSE IF( n.GT.0 ) THEN
723 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
724 $ ten*ulp*temp3, ten*rtunfl )
725 END IF
726 ELSE
727 temp3 = zero
728 vl = zero
729 vu = one
730 END IF
731*
732 CALL zheevx( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
733 $ abstol, m, wa1, z, ldu, work, lwork, rwork,
734 $ iwork, iwork( 5*n+1 ), iinfo )
735 IF( iinfo.NE.0 ) THEN
736 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,A,' // uplo //
737 $ ')', iinfo, n, jtype, ioldsd
738 info = abs( iinfo )
739 IF( iinfo.LT.0 ) THEN
740 RETURN
741 ELSE
742 result( ntest ) = ulpinv
743 result( ntest+1 ) = ulpinv
744 result( ntest+2 ) = ulpinv
745 GO TO 150
746 END IF
747 END IF
748*
749* Do tests 4 and 5.
750*
751 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
752*
753 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
754 $ ldu, tau, work, rwork, result( ntest ) )
755*
756 ntest = ntest + 2
757 CALL zheevx( 'N', 'A', uplo, n, a, ldu, vl, vu, il, iu,
758 $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
759 $ iwork, iwork( 5*n+1 ), iinfo )
760 IF( iinfo.NE.0 ) THEN
761 WRITE( nounit, fmt = 9999 )'ZHEEVX(N,A,' // uplo //
762 $ ')', iinfo, n, jtype, ioldsd
763 info = abs( iinfo )
764 IF( iinfo.LT.0 ) THEN
765 RETURN
766 ELSE
767 result( ntest ) = ulpinv
768 GO TO 150
769 END IF
770 END IF
771*
772* Do test 6.
773*
774 temp1 = zero
775 temp2 = zero
776 DO 140 j = 1, n
777 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
778 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
779 140 CONTINUE
780 result( ntest ) = temp2 / max( unfl,
781 $ ulp*max( temp1, temp2 ) )
782*
783 150 CONTINUE
784 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
785*
786 ntest = ntest + 1
787*
788 CALL zheevx( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
789 $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
790 $ iwork, iwork( 5*n+1 ), iinfo )
791 IF( iinfo.NE.0 ) THEN
792 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,I,' // uplo //
793 $ ')', iinfo, n, jtype, ioldsd
794 info = abs( iinfo )
795 IF( iinfo.LT.0 ) THEN
796 RETURN
797 ELSE
798 result( ntest ) = ulpinv
799 GO TO 160
800 END IF
801 END IF
802*
803* Do tests 7 and 8.
804*
805 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
806*
807 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
808 $ v, ldu, tau, work, rwork, result( ntest ) )
809*
810 ntest = ntest + 2
811*
812 CALL zheevx( 'N', 'I', uplo, n, a, ldu, vl, vu, il, iu,
813 $ abstol, m3, wa3, z, ldu, work, lwork, rwork,
814 $ iwork, iwork( 5*n+1 ), iinfo )
815 IF( iinfo.NE.0 ) THEN
816 WRITE( nounit, fmt = 9999 )'ZHEEVX(N,I,' // uplo //
817 $ ')', iinfo, n, jtype, ioldsd
818 info = abs( iinfo )
819 IF( iinfo.LT.0 ) THEN
820 RETURN
821 ELSE
822 result( ntest ) = ulpinv
823 GO TO 160
824 END IF
825 END IF
826*
827* Do test 9.
828*
829 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
830 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
831 IF( n.GT.0 ) THEN
832 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
833 ELSE
834 temp3 = zero
835 END IF
836 result( ntest ) = ( temp1+temp2 ) /
837 $ max( unfl, temp3*ulp )
838*
839 160 CONTINUE
840 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
841*
842 ntest = ntest + 1
843*
844 CALL zheevx( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
845 $ abstol, m2, wa2, z, ldu, work, lwork, rwork,
846 $ iwork, iwork( 5*n+1 ), iinfo )
847 IF( iinfo.NE.0 ) THEN
848 WRITE( nounit, fmt = 9999 )'ZHEEVX(V,V,' // uplo //
849 $ ')', iinfo, n, jtype, ioldsd
850 info = abs( iinfo )
851 IF( iinfo.LT.0 ) THEN
852 RETURN
853 ELSE
854 result( ntest ) = ulpinv
855 GO TO 170
856 END IF
857 END IF
858*
859* Do tests 10 and 11.
860*
861 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
862*
863 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
864 $ v, ldu, tau, work, rwork, result( ntest ) )
865*
866 ntest = ntest + 2
867*
868 CALL zheevx( 'N', 'V', uplo, n, a, ldu, vl, vu, il, iu,
869 $ abstol, m3, wa3, z, ldu, work, lwork, rwork,
870 $ iwork, iwork( 5*n+1 ), iinfo )
871 IF( iinfo.NE.0 ) THEN
872 WRITE( nounit, fmt = 9999 )'ZHEEVX(N,V,' // uplo //
873 $ ')', iinfo, n, jtype, ioldsd
874 info = abs( iinfo )
875 IF( iinfo.LT.0 ) THEN
876 RETURN
877 ELSE
878 result( ntest ) = ulpinv
879 GO TO 170
880 END IF
881 END IF
882*
883 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
884 result( ntest ) = ulpinv
885 GO TO 170
886 END IF
887*
888* Do test 12.
889*
890 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
891 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
892 IF( n.GT.0 ) THEN
893 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
894 ELSE
895 temp3 = zero
896 END IF
897 result( ntest ) = ( temp1+temp2 ) /
898 $ max( unfl, temp3*ulp )
899*
900 170 CONTINUE
901*
902* Call ZHPEVD and CHPEVX.
903*
904 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
905*
906* Load array WORK with the upper or lower triangular
907* part of the matrix in packed form.
908*
909 IF( iuplo.EQ.1 ) THEN
910 indx = 1
911 DO 190 j = 1, n
912 DO 180 i = 1, j
913 work( indx ) = a( i, j )
914 indx = indx + 1
915 180 CONTINUE
916 190 CONTINUE
917 ELSE
918 indx = 1
919 DO 210 j = 1, n
920 DO 200 i = j, n
921 work( indx ) = a( i, j )
922 indx = indx + 1
923 200 CONTINUE
924 210 CONTINUE
925 END IF
926*
927 ntest = ntest + 1
928 indwrk = n*( n+1 ) / 2 + 1
929 CALL zhpevd( 'V', uplo, n, work, d1, z, ldu,
930 $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
931 $ liwedc, iinfo )
932 IF( iinfo.NE.0 ) THEN
933 WRITE( nounit, fmt = 9999 )'ZHPEVD(V,' // uplo //
934 $ ')', iinfo, n, jtype, ioldsd
935 info = abs( iinfo )
936 IF( iinfo.LT.0 ) THEN
937 RETURN
938 ELSE
939 result( ntest ) = ulpinv
940 result( ntest+1 ) = ulpinv
941 result( ntest+2 ) = ulpinv
942 GO TO 270
943 END IF
944 END IF
945*
946* Do tests 13 and 14.
947*
948 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
949 $ ldu, tau, work, rwork, result( ntest ) )
950*
951 IF( iuplo.EQ.1 ) THEN
952 indx = 1
953 DO 230 j = 1, n
954 DO 220 i = 1, j
955 work( indx ) = a( i, j )
956 indx = indx + 1
957 220 CONTINUE
958 230 CONTINUE
959 ELSE
960 indx = 1
961 DO 250 j = 1, n
962 DO 240 i = j, n
963 work( indx ) = a( i, j )
964 indx = indx + 1
965 240 CONTINUE
966 250 CONTINUE
967 END IF
968*
969 ntest = ntest + 2
970 indwrk = n*( n+1 ) / 2 + 1
971 CALL zhpevd( 'N', uplo, n, work, d3, z, ldu,
972 $ work( indwrk ), lwedc, rwork, lrwedc, iwork,
973 $ liwedc, iinfo )
974 IF( iinfo.NE.0 ) THEN
975 WRITE( nounit, fmt = 9999 )'ZHPEVD(N,' // uplo //
976 $ ')', iinfo, n, jtype, ioldsd
977 info = abs( iinfo )
978 IF( iinfo.LT.0 ) THEN
979 RETURN
980 ELSE
981 result( ntest ) = ulpinv
982 GO TO 270
983 END IF
984 END IF
985*
986* Do test 15.
987*
988 temp1 = zero
989 temp2 = zero
990 DO 260 j = 1, n
991 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
992 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
993 260 CONTINUE
994 result( ntest ) = temp2 / max( unfl,
995 $ ulp*max( temp1, temp2 ) )
996*
997* Load array WORK with the upper or lower triangular part
998* of the matrix in packed form.
999*
1000 270 CONTINUE
1001 IF( iuplo.EQ.1 ) THEN
1002 indx = 1
1003 DO 290 j = 1, n
1004 DO 280 i = 1, j
1005 work( indx ) = a( i, j )
1006 indx = indx + 1
1007 280 CONTINUE
1008 290 CONTINUE
1009 ELSE
1010 indx = 1
1011 DO 310 j = 1, n
1012 DO 300 i = j, n
1013 work( indx ) = a( i, j )
1014 indx = indx + 1
1015 300 CONTINUE
1016 310 CONTINUE
1017 END IF
1018*
1019 ntest = ntest + 1
1020*
1021 IF( n.GT.0 ) THEN
1022 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1023 IF( il.NE.1 ) THEN
1024 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1025 $ ten*ulp*temp3, ten*rtunfl )
1026 ELSE IF( n.GT.0 ) THEN
1027 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1028 $ ten*ulp*temp3, ten*rtunfl )
1029 END IF
1030 IF( iu.NE.n ) THEN
1031 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1032 $ ten*ulp*temp3, ten*rtunfl )
1033 ELSE IF( n.GT.0 ) THEN
1034 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1035 $ ten*ulp*temp3, ten*rtunfl )
1036 END IF
1037 ELSE
1038 temp3 = zero
1039 vl = zero
1040 vu = one
1041 END IF
1042*
1043 CALL zhpevx( 'V', 'A', uplo, n, work, vl, vu, il, iu,
1044 $ abstol, m, wa1, z, ldu, v, rwork, iwork,
1045 $ iwork( 5*n+1 ), iinfo )
1046 IF( iinfo.NE.0 ) THEN
1047 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,A,' // uplo //
1048 $ ')', iinfo, n, jtype, ioldsd
1049 info = abs( iinfo )
1050 IF( iinfo.LT.0 ) THEN
1051 RETURN
1052 ELSE
1053 result( ntest ) = ulpinv
1054 result( ntest+1 ) = ulpinv
1055 result( ntest+2 ) = ulpinv
1056 GO TO 370
1057 END IF
1058 END IF
1059*
1060* Do tests 16 and 17.
1061*
1062 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1063 $ ldu, tau, work, rwork, result( ntest ) )
1064*
1065 ntest = ntest + 2
1066*
1067 IF( iuplo.EQ.1 ) THEN
1068 indx = 1
1069 DO 330 j = 1, n
1070 DO 320 i = 1, j
1071 work( indx ) = a( i, j )
1072 indx = indx + 1
1073 320 CONTINUE
1074 330 CONTINUE
1075 ELSE
1076 indx = 1
1077 DO 350 j = 1, n
1078 DO 340 i = j, n
1079 work( indx ) = a( i, j )
1080 indx = indx + 1
1081 340 CONTINUE
1082 350 CONTINUE
1083 END IF
1084*
1085 CALL zhpevx( 'N', 'A', uplo, n, work, vl, vu, il, iu,
1086 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1087 $ iwork( 5*n+1 ), iinfo )
1088 IF( iinfo.NE.0 ) THEN
1089 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,A,' // uplo //
1090 $ ')', iinfo, n, jtype, ioldsd
1091 info = abs( iinfo )
1092 IF( iinfo.LT.0 ) THEN
1093 RETURN
1094 ELSE
1095 result( ntest ) = ulpinv
1096 GO TO 370
1097 END IF
1098 END IF
1099*
1100* Do test 18.
1101*
1102 temp1 = zero
1103 temp2 = zero
1104 DO 360 j = 1, n
1105 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1106 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1107 360 CONTINUE
1108 result( ntest ) = temp2 / max( unfl,
1109 $ ulp*max( temp1, temp2 ) )
1110*
1111 370 CONTINUE
1112 ntest = ntest + 1
1113 IF( iuplo.EQ.1 ) THEN
1114 indx = 1
1115 DO 390 j = 1, n
1116 DO 380 i = 1, j
1117 work( indx ) = a( i, j )
1118 indx = indx + 1
1119 380 CONTINUE
1120 390 CONTINUE
1121 ELSE
1122 indx = 1
1123 DO 410 j = 1, n
1124 DO 400 i = j, n
1125 work( indx ) = a( i, j )
1126 indx = indx + 1
1127 400 CONTINUE
1128 410 CONTINUE
1129 END IF
1130*
1131 CALL zhpevx( 'V', 'I', uplo, n, work, vl, vu, il, iu,
1132 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1133 $ iwork( 5*n+1 ), iinfo )
1134 IF( iinfo.NE.0 ) THEN
1135 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,I,' // uplo //
1136 $ ')', iinfo, n, jtype, ioldsd
1137 info = abs( iinfo )
1138 IF( iinfo.LT.0 ) THEN
1139 RETURN
1140 ELSE
1141 result( ntest ) = ulpinv
1142 result( ntest+1 ) = ulpinv
1143 result( ntest+2 ) = ulpinv
1144 GO TO 460
1145 END IF
1146 END IF
1147*
1148* Do tests 19 and 20.
1149*
1150 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1151 $ v, ldu, tau, work, rwork, result( ntest ) )
1152*
1153 ntest = ntest + 2
1154*
1155 IF( iuplo.EQ.1 ) THEN
1156 indx = 1
1157 DO 430 j = 1, n
1158 DO 420 i = 1, j
1159 work( indx ) = a( i, j )
1160 indx = indx + 1
1161 420 CONTINUE
1162 430 CONTINUE
1163 ELSE
1164 indx = 1
1165 DO 450 j = 1, n
1166 DO 440 i = j, n
1167 work( indx ) = a( i, j )
1168 indx = indx + 1
1169 440 CONTINUE
1170 450 CONTINUE
1171 END IF
1172*
1173 CALL zhpevx( 'N', 'I', uplo, n, work, vl, vu, il, iu,
1174 $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1175 $ iwork( 5*n+1 ), iinfo )
1176 IF( iinfo.NE.0 ) THEN
1177 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,I,' // uplo //
1178 $ ')', iinfo, n, jtype, ioldsd
1179 info = abs( iinfo )
1180 IF( iinfo.LT.0 ) THEN
1181 RETURN
1182 ELSE
1183 result( ntest ) = ulpinv
1184 GO TO 460
1185 END IF
1186 END IF
1187*
1188* Do test 21.
1189*
1190 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1191 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1192 IF( n.GT.0 ) THEN
1193 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1194 ELSE
1195 temp3 = zero
1196 END IF
1197 result( ntest ) = ( temp1+temp2 ) /
1198 $ max( unfl, temp3*ulp )
1199*
1200 460 CONTINUE
1201 ntest = ntest + 1
1202 IF( iuplo.EQ.1 ) THEN
1203 indx = 1
1204 DO 480 j = 1, n
1205 DO 470 i = 1, j
1206 work( indx ) = a( i, j )
1207 indx = indx + 1
1208 470 CONTINUE
1209 480 CONTINUE
1210 ELSE
1211 indx = 1
1212 DO 500 j = 1, n
1213 DO 490 i = j, n
1214 work( indx ) = a( i, j )
1215 indx = indx + 1
1216 490 CONTINUE
1217 500 CONTINUE
1218 END IF
1219*
1220 CALL zhpevx( 'V', 'V', uplo, n, work, vl, vu, il, iu,
1221 $ abstol, m2, wa2, z, ldu, v, rwork, iwork,
1222 $ iwork( 5*n+1 ), iinfo )
1223 IF( iinfo.NE.0 ) THEN
1224 WRITE( nounit, fmt = 9999 )'ZHPEVX(V,V,' // uplo //
1225 $ ')', iinfo, n, jtype, ioldsd
1226 info = abs( iinfo )
1227 IF( iinfo.LT.0 ) THEN
1228 RETURN
1229 ELSE
1230 result( ntest ) = ulpinv
1231 result( ntest+1 ) = ulpinv
1232 result( ntest+2 ) = ulpinv
1233 GO TO 550
1234 END IF
1235 END IF
1236*
1237* Do tests 22 and 23.
1238*
1239 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1240 $ v, ldu, tau, work, rwork, result( ntest ) )
1241*
1242 ntest = ntest + 2
1243*
1244 IF( iuplo.EQ.1 ) THEN
1245 indx = 1
1246 DO 520 j = 1, n
1247 DO 510 i = 1, j
1248 work( indx ) = a( i, j )
1249 indx = indx + 1
1250 510 CONTINUE
1251 520 CONTINUE
1252 ELSE
1253 indx = 1
1254 DO 540 j = 1, n
1255 DO 530 i = j, n
1256 work( indx ) = a( i, j )
1257 indx = indx + 1
1258 530 CONTINUE
1259 540 CONTINUE
1260 END IF
1261*
1262 CALL zhpevx( 'N', 'V', uplo, n, work, vl, vu, il, iu,
1263 $ abstol, m3, wa3, z, ldu, v, rwork, iwork,
1264 $ iwork( 5*n+1 ), iinfo )
1265 IF( iinfo.NE.0 ) THEN
1266 WRITE( nounit, fmt = 9999 )'ZHPEVX(N,V,' // uplo //
1267 $ ')', iinfo, n, jtype, ioldsd
1268 info = abs( iinfo )
1269 IF( iinfo.LT.0 ) THEN
1270 RETURN
1271 ELSE
1272 result( ntest ) = ulpinv
1273 GO TO 550
1274 END IF
1275 END IF
1276*
1277 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1278 result( ntest ) = ulpinv
1279 GO TO 550
1280 END IF
1281*
1282* Do test 24.
1283*
1284 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1285 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1286 IF( n.GT.0 ) THEN
1287 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1288 ELSE
1289 temp3 = zero
1290 END IF
1291 result( ntest ) = ( temp1+temp2 ) /
1292 $ max( unfl, temp3*ulp )
1293*
1294 550 CONTINUE
1295*
1296* Call ZHBEVD and CHBEVX.
1297*
1298 IF( jtype.LE.7 ) THEN
1299 kd = 0
1300 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1301 kd = max( n-1, 0 )
1302 ELSE
1303 kd = ihbw
1304 END IF
1305*
1306* Load array V with the upper or lower triangular part
1307* of the matrix in band form.
1308*
1309 IF( iuplo.EQ.1 ) THEN
1310 DO 570 j = 1, n
1311 DO 560 i = max( 1, j-kd ), j
1312 v( kd+1+i-j, j ) = a( i, j )
1313 560 CONTINUE
1314 570 CONTINUE
1315 ELSE
1316 DO 590 j = 1, n
1317 DO 580 i = j, min( n, j+kd )
1318 v( 1+i-j, j ) = a( i, j )
1319 580 CONTINUE
1320 590 CONTINUE
1321 END IF
1322*
1323 ntest = ntest + 1
1324 CALL zhbevd( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1325 $ lwedc, rwork, lrwedc, iwork, liwedc, iinfo )
1326 IF( iinfo.NE.0 ) THEN
1327 WRITE( nounit, fmt = 9998 )'ZHBEVD(V,' // uplo //
1328 $ ')', iinfo, n, kd, jtype, ioldsd
1329 info = abs( iinfo )
1330 IF( iinfo.LT.0 ) THEN
1331 RETURN
1332 ELSE
1333 result( ntest ) = ulpinv
1334 result( ntest+1 ) = ulpinv
1335 result( ntest+2 ) = ulpinv
1336 GO TO 650
1337 END IF
1338 END IF
1339*
1340* Do tests 25 and 26.
1341*
1342 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1343 $ ldu, tau, work, rwork, result( ntest ) )
1344*
1345 IF( iuplo.EQ.1 ) THEN
1346 DO 610 j = 1, n
1347 DO 600 i = max( 1, j-kd ), j
1348 v( kd+1+i-j, j ) = a( i, j )
1349 600 CONTINUE
1350 610 CONTINUE
1351 ELSE
1352 DO 630 j = 1, n
1353 DO 620 i = j, min( n, j+kd )
1354 v( 1+i-j, j ) = a( i, j )
1355 620 CONTINUE
1356 630 CONTINUE
1357 END IF
1358*
1359 ntest = ntest + 2
1360 CALL zhbevd( 'N', uplo, n, kd, v, ldu, d3, z, ldu, work,
1361 $ lwedc, rwork, lrwedc, iwork, liwedc, iinfo )
1362 IF( iinfo.NE.0 ) THEN
1363 WRITE( nounit, fmt = 9998 )'ZHBEVD(N,' // uplo //
1364 $ ')', iinfo, n, kd, jtype, ioldsd
1365 info = abs( iinfo )
1366 IF( iinfo.LT.0 ) THEN
1367 RETURN
1368 ELSE
1369 result( ntest ) = ulpinv
1370 GO TO 650
1371 END IF
1372 END IF
1373*
1374* Do test 27.
1375*
1376 temp1 = zero
1377 temp2 = zero
1378 DO 640 j = 1, n
1379 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1380 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1381 640 CONTINUE
1382 result( ntest ) = temp2 / max( unfl,
1383 $ ulp*max( temp1, temp2 ) )
1384*
1385* Load array V with the upper or lower triangular part
1386* of the matrix in band form.
1387*
1388 650 CONTINUE
1389 IF( iuplo.EQ.1 ) THEN
1390 DO 670 j = 1, n
1391 DO 660 i = max( 1, j-kd ), j
1392 v( kd+1+i-j, j ) = a( i, j )
1393 660 CONTINUE
1394 670 CONTINUE
1395 ELSE
1396 DO 690 j = 1, n
1397 DO 680 i = j, min( n, j+kd )
1398 v( 1+i-j, j ) = a( i, j )
1399 680 CONTINUE
1400 690 CONTINUE
1401 END IF
1402*
1403 ntest = ntest + 1
1404 CALL zhbevx( 'V', 'A', uplo, n, kd, v, ldu, u, ldu, vl,
1405 $ vu, il, iu, abstol, m, wa1, z, ldu, work,
1406 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1407 IF( iinfo.NE.0 ) THEN
1408 WRITE( nounit, fmt = 9999 )'ZHBEVX(V,A,' // uplo //
1409 $ ')', iinfo, n, kd, jtype, ioldsd
1410 info = abs( iinfo )
1411 IF( iinfo.LT.0 ) THEN
1412 RETURN
1413 ELSE
1414 result( ntest ) = ulpinv
1415 result( ntest+1 ) = ulpinv
1416 result( ntest+2 ) = ulpinv
1417 GO TO 750
1418 END IF
1419 END IF
1420*
1421* Do tests 28 and 29.
1422*
1423 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1424 $ ldu, tau, work, rwork, result( ntest ) )
1425*
1426 ntest = ntest + 2
1427*
1428 IF( iuplo.EQ.1 ) THEN
1429 DO 710 j = 1, n
1430 DO 700 i = max( 1, j-kd ), j
1431 v( kd+1+i-j, j ) = a( i, j )
1432 700 CONTINUE
1433 710 CONTINUE
1434 ELSE
1435 DO 730 j = 1, n
1436 DO 720 i = j, min( n, j+kd )
1437 v( 1+i-j, j ) = a( i, j )
1438 720 CONTINUE
1439 730 CONTINUE
1440 END IF
1441*
1442 CALL zhbevx( 'N', 'A', uplo, n, kd, v, ldu, u, ldu, vl,
1443 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1444 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1445 IF( iinfo.NE.0 ) THEN
1446 WRITE( nounit, fmt = 9998 )'ZHBEVX(N,A,' // uplo //
1447 $ ')', iinfo, n, kd, jtype, ioldsd
1448 info = abs( iinfo )
1449 IF( iinfo.LT.0 ) THEN
1450 RETURN
1451 ELSE
1452 result( ntest ) = ulpinv
1453 GO TO 750
1454 END IF
1455 END IF
1456*
1457* Do test 30.
1458*
1459 temp1 = zero
1460 temp2 = zero
1461 DO 740 j = 1, n
1462 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1463 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1464 740 CONTINUE
1465 result( ntest ) = temp2 / max( unfl,
1466 $ ulp*max( temp1, temp2 ) )
1467*
1468* Load array V with the upper or lower triangular part
1469* of the matrix in band form.
1470*
1471 750 CONTINUE
1472 ntest = ntest + 1
1473 IF( iuplo.EQ.1 ) THEN
1474 DO 770 j = 1, n
1475 DO 760 i = max( 1, j-kd ), j
1476 v( kd+1+i-j, j ) = a( i, j )
1477 760 CONTINUE
1478 770 CONTINUE
1479 ELSE
1480 DO 790 j = 1, n
1481 DO 780 i = j, min( n, j+kd )
1482 v( 1+i-j, j ) = a( i, j )
1483 780 CONTINUE
1484 790 CONTINUE
1485 END IF
1486*
1487 CALL zhbevx( 'V', 'I', uplo, n, kd, v, ldu, u, ldu, vl,
1488 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1489 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1490 IF( iinfo.NE.0 ) THEN
1491 WRITE( nounit, fmt = 9998 )'ZHBEVX(V,I,' // uplo //
1492 $ ')', iinfo, n, kd, jtype, ioldsd
1493 info = abs( iinfo )
1494 IF( iinfo.LT.0 ) THEN
1495 RETURN
1496 ELSE
1497 result( ntest ) = ulpinv
1498 result( ntest+1 ) = ulpinv
1499 result( ntest+2 ) = ulpinv
1500 GO TO 840
1501 END IF
1502 END IF
1503*
1504* Do tests 31 and 32.
1505*
1506 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1507 $ v, ldu, tau, work, rwork, result( ntest ) )
1508*
1509 ntest = ntest + 2
1510*
1511 IF( iuplo.EQ.1 ) THEN
1512 DO 810 j = 1, n
1513 DO 800 i = max( 1, j-kd ), j
1514 v( kd+1+i-j, j ) = a( i, j )
1515 800 CONTINUE
1516 810 CONTINUE
1517 ELSE
1518 DO 830 j = 1, n
1519 DO 820 i = j, min( n, j+kd )
1520 v( 1+i-j, j ) = a( i, j )
1521 820 CONTINUE
1522 830 CONTINUE
1523 END IF
1524 CALL zhbevx( 'N', 'I', uplo, n, kd, v, ldu, u, ldu, vl,
1525 $ vu, il, iu, abstol, m3, wa3, z, ldu, work,
1526 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1527 IF( iinfo.NE.0 ) THEN
1528 WRITE( nounit, fmt = 9998 )'ZHBEVX(N,I,' // uplo //
1529 $ ')', iinfo, n, kd, jtype, ioldsd
1530 info = abs( iinfo )
1531 IF( iinfo.LT.0 ) THEN
1532 RETURN
1533 ELSE
1534 result( ntest ) = ulpinv
1535 GO TO 840
1536 END IF
1537 END IF
1538*
1539* Do test 33.
1540*
1541 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1542 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1543 IF( n.GT.0 ) THEN
1544 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1545 ELSE
1546 temp3 = zero
1547 END IF
1548 result( ntest ) = ( temp1+temp2 ) /
1549 $ max( unfl, temp3*ulp )
1550*
1551* Load array V with the upper or lower triangular part
1552* of the matrix in band form.
1553*
1554 840 CONTINUE
1555 ntest = ntest + 1
1556 IF( iuplo.EQ.1 ) THEN
1557 DO 860 j = 1, n
1558 DO 850 i = max( 1, j-kd ), j
1559 v( kd+1+i-j, j ) = a( i, j )
1560 850 CONTINUE
1561 860 CONTINUE
1562 ELSE
1563 DO 880 j = 1, n
1564 DO 870 i = j, min( n, j+kd )
1565 v( 1+i-j, j ) = a( i, j )
1566 870 CONTINUE
1567 880 CONTINUE
1568 END IF
1569 CALL zhbevx( 'V', 'V', uplo, n, kd, v, ldu, u, ldu, vl,
1570 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
1571 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1572 IF( iinfo.NE.0 ) THEN
1573 WRITE( nounit, fmt = 9998 )'ZHBEVX(V,V,' // uplo //
1574 $ ')', iinfo, n, kd, jtype, ioldsd
1575 info = abs( iinfo )
1576 IF( iinfo.LT.0 ) THEN
1577 RETURN
1578 ELSE
1579 result( ntest ) = ulpinv
1580 result( ntest+1 ) = ulpinv
1581 result( ntest+2 ) = ulpinv
1582 GO TO 930
1583 END IF
1584 END IF
1585*
1586* Do tests 34 and 35.
1587*
1588 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1589 $ v, ldu, tau, work, rwork, result( ntest ) )
1590*
1591 ntest = ntest + 2
1592*
1593 IF( iuplo.EQ.1 ) THEN
1594 DO 900 j = 1, n
1595 DO 890 i = max( 1, j-kd ), j
1596 v( kd+1+i-j, j ) = a( i, j )
1597 890 CONTINUE
1598 900 CONTINUE
1599 ELSE
1600 DO 920 j = 1, n
1601 DO 910 i = j, min( n, j+kd )
1602 v( 1+i-j, j ) = a( i, j )
1603 910 CONTINUE
1604 920 CONTINUE
1605 END IF
1606 CALL zhbevx( 'N', 'V', uplo, n, kd, v, ldu, u, ldu, vl,
1607 $ vu, il, iu, abstol, m3, wa3, z, ldu, work,
1608 $ rwork, iwork, iwork( 5*n+1 ), iinfo )
1609 IF( iinfo.NE.0 ) THEN
1610 WRITE( nounit, fmt = 9998 )'ZHBEVX(N,V,' // uplo //
1611 $ ')', iinfo, n, kd, jtype, ioldsd
1612 info = abs( iinfo )
1613 IF( iinfo.LT.0 ) THEN
1614 RETURN
1615 ELSE
1616 result( ntest ) = ulpinv
1617 GO TO 930
1618 END IF
1619 END IF
1620*
1621 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1622 result( ntest ) = ulpinv
1623 GO TO 930
1624 END IF
1625*
1626* Do test 36.
1627*
1628 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1629 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1630 IF( n.GT.0 ) THEN
1631 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1632 ELSE
1633 temp3 = zero
1634 END IF
1635 result( ntest ) = ( temp1+temp2 ) /
1636 $ max( unfl, temp3*ulp )
1637*
1638 930 CONTINUE
1639*
1640* Call ZHEEV
1641*
1642 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
1643*
1644 ntest = ntest + 1
1645 CALL zheev( 'V', uplo, n, a, ldu, d1, work, lwork, rwork,
1646 $ iinfo )
1647 IF( iinfo.NE.0 ) THEN
1648 WRITE( nounit, fmt = 9999 )'ZHEEV(V,' // uplo // ')',
1649 $ iinfo, n, jtype, ioldsd
1650 info = abs( iinfo )
1651 IF( iinfo.LT.0 ) THEN
1652 RETURN
1653 ELSE
1654 result( ntest ) = ulpinv
1655 result( ntest+1 ) = ulpinv
1656 result( ntest+2 ) = ulpinv
1657 GO TO 950
1658 END IF
1659 END IF
1660*
1661* Do tests 37 and 38
1662*
1663 CALL zhet21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
1664 $ ldu, tau, work, rwork, result( ntest ) )
1665*
1666 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1667*
1668 ntest = ntest + 2
1669 CALL zheev( 'N', uplo, n, a, ldu, d3, work, lwork, rwork,
1670 $ iinfo )
1671 IF( iinfo.NE.0 ) THEN
1672 WRITE( nounit, fmt = 9999 )'ZHEEV(N,' // uplo // ')',
1673 $ iinfo, n, jtype, ioldsd
1674 info = abs( iinfo )
1675 IF( iinfo.LT.0 ) THEN
1676 RETURN
1677 ELSE
1678 result( ntest ) = ulpinv
1679 GO TO 950
1680 END IF
1681 END IF
1682*
1683* Do test 39
1684*
1685 temp1 = zero
1686 temp2 = zero
1687 DO 940 j = 1, n
1688 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1689 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1690 940 CONTINUE
1691 result( ntest ) = temp2 / max( unfl,
1692 $ ulp*max( temp1, temp2 ) )
1693*
1694 950 CONTINUE
1695*
1696 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1697*
1698* Call ZHPEV
1699*
1700* Load array WORK with the upper or lower triangular
1701* part of the matrix in packed form.
1702*
1703 IF( iuplo.EQ.1 ) THEN
1704 indx = 1
1705 DO 970 j = 1, n
1706 DO 960 i = 1, j
1707 work( indx ) = a( i, j )
1708 indx = indx + 1
1709 960 CONTINUE
1710 970 CONTINUE
1711 ELSE
1712 indx = 1
1713 DO 990 j = 1, n
1714 DO 980 i = j, n
1715 work( indx ) = a( i, j )
1716 indx = indx + 1
1717 980 CONTINUE
1718 990 CONTINUE
1719 END IF
1720*
1721 ntest = ntest + 1
1722 indwrk = n*( n+1 ) / 2 + 1
1723 CALL zhpev( 'V', uplo, n, work, d1, z, ldu,
1724 $ work( indwrk ), rwork, iinfo )
1725 IF( iinfo.NE.0 ) THEN
1726 WRITE( nounit, fmt = 9999 )'ZHPEV(V,' // uplo // ')',
1727 $ iinfo, n, jtype, ioldsd
1728 info = abs( iinfo )
1729 IF( iinfo.LT.0 ) THEN
1730 RETURN
1731 ELSE
1732 result( ntest ) = ulpinv
1733 result( ntest+1 ) = ulpinv
1734 result( ntest+2 ) = ulpinv
1735 GO TO 1050
1736 END IF
1737 END IF
1738*
1739* Do tests 40 and 41.
1740*
1741 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1742 $ ldu, tau, work, rwork, result( ntest ) )
1743*
1744 IF( iuplo.EQ.1 ) THEN
1745 indx = 1
1746 DO 1010 j = 1, n
1747 DO 1000 i = 1, j
1748 work( indx ) = a( i, j )
1749 indx = indx + 1
1750 1000 CONTINUE
1751 1010 CONTINUE
1752 ELSE
1753 indx = 1
1754 DO 1030 j = 1, n
1755 DO 1020 i = j, n
1756 work( indx ) = a( i, j )
1757 indx = indx + 1
1758 1020 CONTINUE
1759 1030 CONTINUE
1760 END IF
1761*
1762 ntest = ntest + 2
1763 indwrk = n*( n+1 ) / 2 + 1
1764 CALL zhpev( 'N', uplo, n, work, d3, z, ldu,
1765 $ work( indwrk ), rwork, iinfo )
1766 IF( iinfo.NE.0 ) THEN
1767 WRITE( nounit, fmt = 9999 )'ZHPEV(N,' // uplo // ')',
1768 $ iinfo, n, jtype, ioldsd
1769 info = abs( iinfo )
1770 IF( iinfo.LT.0 ) THEN
1771 RETURN
1772 ELSE
1773 result( ntest ) = ulpinv
1774 GO TO 1050
1775 END IF
1776 END IF
1777*
1778* Do test 42
1779*
1780 temp1 = zero
1781 temp2 = zero
1782 DO 1040 j = 1, n
1783 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1784 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1785 1040 CONTINUE
1786 result( ntest ) = temp2 / max( unfl,
1787 $ ulp*max( temp1, temp2 ) )
1788*
1789 1050 CONTINUE
1790*
1791* Call ZHBEV
1792*
1793 IF( jtype.LE.7 ) THEN
1794 kd = 0
1795 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
1796 kd = max( n-1, 0 )
1797 ELSE
1798 kd = ihbw
1799 END IF
1800*
1801* Load array V with the upper or lower triangular part
1802* of the matrix in band form.
1803*
1804 IF( iuplo.EQ.1 ) THEN
1805 DO 1070 j = 1, n
1806 DO 1060 i = max( 1, j-kd ), j
1807 v( kd+1+i-j, j ) = a( i, j )
1808 1060 CONTINUE
1809 1070 CONTINUE
1810 ELSE
1811 DO 1090 j = 1, n
1812 DO 1080 i = j, min( n, j+kd )
1813 v( 1+i-j, j ) = a( i, j )
1814 1080 CONTINUE
1815 1090 CONTINUE
1816 END IF
1817*
1818 ntest = ntest + 1
1819 CALL zhbev( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
1820 $ rwork, iinfo )
1821 IF( iinfo.NE.0 ) THEN
1822 WRITE( nounit, fmt = 9998 )'ZHBEV(V,' // uplo // ')',
1823 $ iinfo, n, kd, jtype, ioldsd
1824 info = abs( iinfo )
1825 IF( iinfo.LT.0 ) THEN
1826 RETURN
1827 ELSE
1828 result( ntest ) = ulpinv
1829 result( ntest+1 ) = ulpinv
1830 result( ntest+2 ) = ulpinv
1831 GO TO 1140
1832 END IF
1833 END IF
1834*
1835* Do tests 43 and 44.
1836*
1837 CALL zhet21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1838 $ ldu, tau, work, rwork, result( ntest ) )
1839*
1840 IF( iuplo.EQ.1 ) THEN
1841 DO 1110 j = 1, n
1842 DO 1100 i = max( 1, j-kd ), j
1843 v( kd+1+i-j, j ) = a( i, j )
1844 1100 CONTINUE
1845 1110 CONTINUE
1846 ELSE
1847 DO 1130 j = 1, n
1848 DO 1120 i = j, min( n, j+kd )
1849 v( 1+i-j, j ) = a( i, j )
1850 1120 CONTINUE
1851 1130 CONTINUE
1852 END IF
1853*
1854 ntest = ntest + 2
1855 CALL zhbev( 'N', uplo, n, kd, v, ldu, d3, z, ldu, work,
1856 $ rwork, iinfo )
1857 IF( iinfo.NE.0 ) THEN
1858 WRITE( nounit, fmt = 9998 )'ZHBEV(N,' // uplo // ')',
1859 $ iinfo, n, kd, jtype, ioldsd
1860 info = abs( iinfo )
1861 IF( iinfo.LT.0 ) THEN
1862 RETURN
1863 ELSE
1864 result( ntest ) = ulpinv
1865 GO TO 1140
1866 END IF
1867 END IF
1868*
1869 1140 CONTINUE
1870*
1871* Do test 45.
1872*
1873 temp1 = zero
1874 temp2 = zero
1875 DO 1150 j = 1, n
1876 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1877 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1878 1150 CONTINUE
1879 result( ntest ) = temp2 / max( unfl,
1880 $ ulp*max( temp1, temp2 ) )
1881*
1882 CALL zlacpy( ' ', n, n, a, lda, v, ldu )
1883 ntest = ntest + 1
1884 CALL zheevr( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
1885 $ abstol, m, wa1, z, ldu, iwork, work, lwork,
1886 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1887 $ iinfo )
1888 IF( iinfo.NE.0 ) THEN
1889 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,A,' // uplo //
1890 $ ')', iinfo, n, jtype, ioldsd
1891 info = abs( iinfo )
1892 IF( iinfo.LT.0 ) THEN
1893 RETURN
1894 ELSE
1895 result( ntest ) = ulpinv
1896 result( ntest+1 ) = ulpinv
1897 result( ntest+2 ) = ulpinv
1898 GO TO 1170
1899 END IF
1900 END IF
1901*
1902* Do tests 45 and 46 (or ... )
1903*
1904 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1905*
1906 CALL zhet21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1907 $ ldu, tau, work, rwork, result( ntest ) )
1908*
1909 ntest = ntest + 2
1910 CALL zheevr( 'N', 'A', uplo, n, a, ldu, vl, vu, il, iu,
1911 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
1912 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1913 $ iinfo )
1914 IF( iinfo.NE.0 ) THEN
1915 WRITE( nounit, fmt = 9999 )'ZHEEVR(N,A,' // uplo //
1916 $ ')', iinfo, n, jtype, ioldsd
1917 info = abs( iinfo )
1918 IF( iinfo.LT.0 ) THEN
1919 RETURN
1920 ELSE
1921 result( ntest ) = ulpinv
1922 GO TO 1170
1923 END IF
1924 END IF
1925*
1926* Do test 47 (or ... )
1927*
1928 temp1 = zero
1929 temp2 = zero
1930 DO 1160 j = 1, n
1931 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1932 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1933 1160 CONTINUE
1934 result( ntest ) = temp2 / max( unfl,
1935 $ ulp*max( temp1, temp2 ) )
1936*
1937 1170 CONTINUE
1938*
1939 ntest = ntest + 1
1940 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1941 CALL zheevr( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
1942 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
1943 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1944 $ iinfo )
1945 IF( iinfo.NE.0 ) THEN
1946 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,I,' // uplo //
1947 $ ')', iinfo, n, jtype, ioldsd
1948 info = abs( iinfo )
1949 IF( iinfo.LT.0 ) THEN
1950 RETURN
1951 ELSE
1952 result( ntest ) = ulpinv
1953 result( ntest+1 ) = ulpinv
1954 result( ntest+2 ) = ulpinv
1955 GO TO 1180
1956 END IF
1957 END IF
1958*
1959* Do tests 48 and 49 (or +??)
1960*
1961 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1962*
1963 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1964 $ v, ldu, tau, work, rwork, result( ntest ) )
1965*
1966 ntest = ntest + 2
1967 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1968 CALL zheevr( 'N', 'I', uplo, n, a, ldu, vl, vu, il, iu,
1969 $ abstol, m3, wa3, z, ldu, iwork, work, lwork,
1970 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1971 $ iinfo )
1972 IF( iinfo.NE.0 ) THEN
1973 WRITE( nounit, fmt = 9999 )'ZHEEVR(N,I,' // uplo //
1974 $ ')', iinfo, n, jtype, ioldsd
1975 info = abs( iinfo )
1976 IF( iinfo.LT.0 ) THEN
1977 RETURN
1978 ELSE
1979 result( ntest ) = ulpinv
1980 GO TO 1180
1981 END IF
1982 END IF
1983*
1984* Do test 50 (or +??)
1985*
1986 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1987 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1988 result( ntest ) = ( temp1+temp2 ) /
1989 $ max( unfl, ulp*temp3 )
1990 1180 CONTINUE
1991*
1992 ntest = ntest + 1
1993 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
1994 CALL zheevr( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
1995 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
1996 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
1997 $ iinfo )
1998 IF( iinfo.NE.0 ) THEN
1999 WRITE( nounit, fmt = 9999 )'ZHEEVR(V,V,' // uplo //
2000 $ ')', iinfo, n, jtype, ioldsd
2001 info = abs( iinfo )
2002 IF( iinfo.LT.0 ) THEN
2003 RETURN
2004 ELSE
2005 result( ntest ) = ulpinv
2006 result( ntest+1 ) = ulpinv
2007 result( ntest+2 ) = ulpinv
2008 GO TO 1190
2009 END IF
2010 END IF
2011*
2012* Do tests 51 and 52 (or +??)
2013*
2014 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2015*
2016 CALL zhet22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2017 $ v, ldu, tau, work, rwork, result( ntest ) )
2018*
2019 ntest = ntest + 2
2020 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2021 CALL zheevr( 'N', 'V', uplo, n, a, ldu, vl, vu, il, iu,
2022 $ abstol, m3, wa3, z, ldu, iwork, work, lwork,
2023 $ rwork, lrwork, iwork( 2*n+1 ), liwork-2*n,
2024 $ iinfo )
2025 IF( iinfo.NE.0 ) THEN
2026 WRITE( nounit, fmt = 9999 )'ZHEEVR(N,V,' // uplo //
2027 $ ')', iinfo, n, jtype, ioldsd
2028 info = abs( iinfo )
2029 IF( iinfo.LT.0 ) THEN
2030 RETURN
2031 ELSE
2032 result( ntest ) = ulpinv
2033 GO TO 1190
2034 END IF
2035 END IF
2036*
2037 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2038 result( ntest ) = ulpinv
2039 GO TO 1190
2040 END IF
2041*
2042* Do test 52 (or +??)
2043*
2044 temp1 = dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2045 temp2 = dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2046 IF( n.GT.0 ) THEN
2047 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2048 ELSE
2049 temp3 = zero
2050 END IF
2051 result( ntest ) = ( temp1+temp2 ) /
2052 $ max( unfl, temp3*ulp )
2053*
2054 CALL zlacpy( ' ', n, n, v, ldu, a, lda )
2055*
2056*
2057*
2058*
2059* Load array V with the upper or lower triangular part
2060* of the matrix in band form.
2061*
2062 1190 CONTINUE
2063*
2064 1200 CONTINUE
2065*
2066* End of Loop -- Check for RESULT(j) > THRESH
2067*
2068 ntestt = ntestt + ntest
2069 CALL dlafts( 'ZST', n, n, jtype, ntest, result, ioldsd,
2070 $ thresh, nounit, nerrs )
2071*
2072 1210 CONTINUE
2073 1220 CONTINUE
2074*
2075* Summary
2076*
2077 CALL alasvm( 'ZST', nounit, nerrs, ntestt, 0 )
2078*
2079 9999 FORMAT( ' ZDRVST: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2080 $ ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
2081 9998 FORMAT( ' ZDRVST: ', a, ' returned INFO=', i6, / 9x, 'N=', i6,
2082 $ ', KD=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
2083 $ ')' )
2084*
2085 RETURN
2086*
2087* End of ZDRVST
2088*
2089 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(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(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(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(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(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(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(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 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 zdrvst(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)
ZDRVST
Definition zdrvst.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