LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zdrvst()

subroutine zdrvst ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
integer  NOUNIT,
complex*16, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  D1,
double precision, dimension( * )  D2,
double precision, dimension( * )  D3,
double precision, dimension( * )  WA1,
double precision, dimension( * )  WA2,
double precision, dimension( * )  WA3,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( ldu, * )  V,
complex*16, dimension( * )  TAU,
complex*16, dimension( ldu, * )  Z,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

ZDRVST

Purpose:
      ZDRVST  checks the Hermitian eigenvalue problem drivers.

              ZHEEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix,
              using a divide-and-conquer algorithm.

              ZHEEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix.

              ZHEEVR computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix
              using the Relatively Robust Representation where it can.

              ZHPEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage, using a divide-and-conquer algorithm.

              ZHPEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage.

              ZHBEVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix,
              using a divide-and-conquer algorithm.

              ZHBEVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix.

              ZHEEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix.

              ZHPEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian matrix in packed
              storage.

              ZHBEV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian band matrix.

      When ZDRVST is called, a number of matrix "sizes" ("n's") and a
      number of matrix "types" are specified.  For each size ("n")
      and each type of matrix, one matrix will be generated and used
      to test the appropriate drivers.  For each matrix and each
      driver routine called, the following tests will be performed:

      (1)     | A - Z D Z' | / ( |A| n ulp )

      (2)     | I - Z Z' | / ( n ulp )

      (3)     | D1 - D2 | / ( |D1| ulp )

      where Z is the matrix of eigenvectors returned when the
      eigenvector option is given and D1 and D2 are the eigenvalues
      returned with and without the eigenvector option.

      The "sizes" are specified by an array NN(1:NSIZES); the value of
      each element NN(j) specifies one size.
      The "types" are specified by a logical array DOTYPE( 1:NTYPES );
      if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
      Currently, the list of possible types is:

      (1)  The zero matrix.
      (2)  The identity matrix.

      (3)  A diagonal matrix with evenly spaced entries
           1, ..., ULP  and random signs.
           (ULP = (first number larger than 1) - 1 )
      (4)  A diagonal matrix with geometrically spaced entries
           1, ..., ULP  and random signs.
      (5)  A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
           and random signs.

      (6)  Same as (4), but multiplied by SQRT( overflow threshold )
      (7)  Same as (4), but multiplied by SQRT( underflow threshold )

      (8)  A matrix of the form  U* D U, where U is unitary and
           D has evenly spaced entries 1, ..., ULP with random signs
           on the diagonal.

      (9)  A matrix of the form  U* D U, where U is unitary and
           D has geometrically spaced entries 1, ..., ULP with random
           signs on the diagonal.

      (10) A matrix of the form  U* D U, where U is unitary and
           D has "clustered" entries 1, ULP,..., ULP with random
           signs on the diagonal.

      (11) Same as (8), but multiplied by SQRT( overflow threshold )
      (12) Same as (8), but multiplied by SQRT( underflow threshold )

      (13) Symmetric matrix with random entries chosen from (-1,1).
      (14) Same as (13), but multiplied by SQRT( overflow threshold )
      (15) Same as (13), but multiplied by SQRT( underflow threshold )
      (16) A band matrix with half bandwidth randomly chosen between
           0 and N-1, with evenly spaced eigenvalues 1, ..., ULP
           with random signs.
      (17) Same as (16), but multiplied by SQRT( overflow threshold )
      (18) Same as (16), but multiplied by SQRT( underflow threshold )
  NSIZES  INTEGER
          The number of sizes of matrices to use.  If it is zero,
          ZDRVST does nothing.  It must be at least zero.
          Not modified.

  NN      INTEGER array, dimension (NSIZES)
          An array containing the sizes to be used for the matrices.
          Zero values will be skipped.  The values must be at least
          zero.
          Not modified.

  NTYPES  INTEGER
          The number of elements in DOTYPE.   If it is zero, ZDRVST
          does nothing.  It must be at least zero.  If it is MAXTYP+1
          and NSIZES is 1, then an additional type, MAXTYP+1 is
          defined, which is to use whatever matrix is in A.  This
          is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
          DOTYPE(MAXTYP+1) is .TRUE. .
          Not modified.

  DOTYPE  LOGICAL array, dimension (NTYPES)
          If DOTYPE(j) is .TRUE., then for each size in NN a
          matrix of that size and of type j will be generated.
          If NTYPES is smaller than the maximum number of types
          defined (PARAMETER MAXTYP), then types NTYPES+1 through
          MAXTYP will not be generated.  If NTYPES is larger
          than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
          will be ignored.
          Not modified.

  ISEED   INTEGER array, dimension (4)
          On entry ISEED specifies the seed of the random number
          generator. The array elements should be between 0 and 4095;
          if not they will be reduced mod 4096.  Also, ISEED(4) must
          be odd.  The random number generator uses a linear
          congruential sequence limited to small integers, and so
          should produce machine independent random numbers. The
          values of ISEED are changed on exit, and can be used in the
          next call to ZDRVST to continue the same random number
          sequence.
          Modified.

  THRESH  DOUBLE PRECISION
          A test will count as "failed" if the "error", computed as
          described above, exceeds THRESH.  Note that the error
          is scaled to be O(1), so THRESH should be a reasonably
          small multiple of 1, e.g., 10 or 100.  In particular,
          it should not depend on the precision (single vs. double)
          or the size of the matrix.  It must be at least zero.
          Not modified.

  NOUNIT  INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
          Not modified.

  A       COMPLEX*16 array, dimension (LDA , max(NN))
          Used to hold the matrix whose eigenvalues are to be
          computed.  On exit, A contains the last matrix actually
          used.
          Modified.

  LDA     INTEGER
          The leading dimension of A.  It must be at
          least 1 and at least max( NN ).
          Not modified.

  D1      DOUBLE PRECISION array, dimension (max(NN))
          The eigenvalues of A, as computed by ZSTEQR simlutaneously
          with Z.  On exit, the eigenvalues in D1 correspond with the
          matrix in A.
          Modified.

  D2      DOUBLE PRECISION array, dimension (max(NN))
          The eigenvalues of A, as computed by ZSTEQR if Z is not
          computed.  On exit, the eigenvalues in D2 correspond with
          the matrix in A.
          Modified.

  D3      DOUBLE PRECISION array, dimension (max(NN))
          The eigenvalues of A, as computed by DSTERF.  On exit, the
          eigenvalues in D3 correspond with the matrix in A.
          Modified.

  WA1     DOUBLE PRECISION array, dimension

  WA2     DOUBLE PRECISION array, dimension

  WA3     DOUBLE PRECISION array, dimension

  U       COMPLEX*16 array, dimension (LDU, max(NN))
          The unitary matrix computed by ZHETRD + ZUNGC3.
          Modified.

  LDU     INTEGER
          The leading dimension of U, Z, and V.  It must be at
          least 1 and at least max( NN ).
          Not modified.

  V       COMPLEX*16 array, dimension (LDU, max(NN))
          The Housholder vectors computed by ZHETRD in reducing A to
          tridiagonal form.
          Modified.

  TAU     COMPLEX*16 array, dimension (max(NN))
          The Householder factors computed by ZHETRD in reducing A
          to tridiagonal form.
          Modified.

  Z       COMPLEX*16 array, dimension (LDU, max(NN))
          The unitary matrix of eigenvectors computed by ZHEEVD,
          ZHEEVX, ZHPEVD, CHPEVX, ZHBEVD, and CHBEVX.
          Modified.

  WORK  - COMPLEX*16 array of dimension ( LWORK )
           Workspace.
           Modified.

  LWORK - INTEGER
           The number of entries in WORK.  This must be at least
           2*max( NN(j), 2 )**2.
           Not modified.

  RWORK   DOUBLE PRECISION array, dimension (3*max(NN))
           Workspace.
           Modified.

  LRWORK - INTEGER
           The number of entries in RWORK.

  IWORK   INTEGER array, dimension (6*max(NN))
          Workspace.
          Modified.

  LIWORK - INTEGER
           The number of entries in IWORK.

  RESULT  DOUBLE PRECISION array, dimension (??)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid
          overflow.
          Modified.

  INFO    INTEGER
          If 0, then everything ran OK.
           -1: NSIZES < 0
           -2: Some NN(j) < 0
           -3: NTYPES < 0
           -5: THRESH < 0
           -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
          -16: LDU < 1 or LDU < NMAX.
          -21: LWORK too small.
          If  DLATMR, SLATMS, ZHETRD, DORGC3, ZSTEQR, DSTERF,
              or DORMC2 returns an error code, the
              absolute value of it is returned.
          Modified.

-----------------------------------------------------------------------

       Some Local Variables and Parameters:
       ---- ----- --------- --- ----------
       ZERO, ONE       Real 0 and 1.
       MAXTYP          The number of types defined.
       NTEST           The number of tests performed, or which can
                       be performed so far, for the current matrix.
       NTESTT          The total number of tests performed so far.
       NMAX            Largest value in NN.
       NMATS           The number of matrices generated so far.
       NERRS           The number of tests which have exceeded THRESH
                       so far (computed by DLAFTS).
       COND, IMODE     Values to be passed to the matrix generators.
       ANORM           Norm of A; passed to matrix generators.

       OVFL, UNFL      Overflow and underflow thresholds.
       ULP, ULPINV     Finest relative precision and its inverse.
       RTOVFL, RTUNFL  Square roots of the previous 2 values.
               The following four arrays decode JTYPE:
       KTYPE(j)        The general type (1-10) for type "j".
       KMODE(j)        The MODE value to be passed to the matrix
                       generator for type "j".
       KMAGN(j)        The order of magnitude ( O(1),
                       O(overflow^(1/2) ), O(underflow^(1/2) )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 334 of file zdrvst.f.

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