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

◆ zdrvsg2stg()

subroutine zdrvsg2stg ( 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,
complex*16, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  D,
double precision, dimension( * )  D2,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( lda, * )  AB,
complex*16, dimension( ldb, * )  BB,
complex*16, dimension( * )  AP,
complex*16, dimension( * )  BP,
complex*16, dimension( * )  WORK,
integer  NWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

ZDRVSG2STG

Purpose:
      ZDRVSG2STG checks the complex Hermitian generalized eigenproblem
      drivers.

              ZHEGV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem.

              ZHEGVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem using a divide and conquer algorithm.

              ZHEGVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem.

              ZHPGV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem in packed storage.

              ZHPGVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem in packed storage using a divide and
              conquer algorithm.

              ZHPGVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite generalized
              eigenproblem in packed storage.

              ZHBGV computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite banded
              generalized eigenproblem.

              ZHBGVD computes all eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite banded
              generalized eigenproblem using a divide and conquer
              algorithm.

              ZHBGVX computes selected eigenvalues and, optionally,
              eigenvectors of a complex Hermitian-definite banded
              generalized eigenproblem.

      When ZDRVSG2STG 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 A of the given type will be
      generated; a random well-conditioned matrix B is also generated
      and the pair (A,B) is used to test the drivers.

      For each pair (A,B), the following tests are performed:

      (1) ZHEGV with ITYPE = 1 and UPLO ='U':

              | A Z - B Z D | / ( |A| |Z| n ulp )
              | D - D2 | / ( |D| ulp )   where D is computed by
                                         ZHEGV and  D2 is computed by
                                         ZHEGV_2STAGE. This test is
                                         only performed for DSYGV

      (2) as (1) but calling ZHPGV
      (3) as (1) but calling ZHBGV
      (4) as (1) but with UPLO = 'L'
      (5) as (4) but calling ZHPGV
      (6) as (4) but calling ZHBGV

      (7) ZHEGV with ITYPE = 2 and UPLO ='U':

              | A B Z - Z D | / ( |A| |Z| n ulp )

      (8) as (7) but calling ZHPGV
      (9) as (7) but with UPLO = 'L'
      (10) as (9) but calling ZHPGV

      (11) ZHEGV with ITYPE = 3 and UPLO ='U':

              | B A Z - Z D | / ( |A| |Z| n ulp )

      (12) as (11) but calling ZHPGV
      (13) as (11) but with UPLO = 'L'
      (14) as (13) but calling ZHPGV

      ZHEGVD, ZHPGVD and ZHBGVD performed the same 14 tests.

      ZHEGVX, ZHPGVX and ZHBGVX performed the above 14 tests with
      the parameter RANGE = 'A', 'N' and 'I', respectively.

      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.
      This type is used for the matrix A which has half-bandwidth KA.
      B is generated as a well-conditioned positive definite matrix
      with half-bandwidth KB (<= KA).
      Currently, the list of possible types for A 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) Hermitian 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) Same as (8), but with KA = 1 and KB = 1
      (17) Same as (8), but with KA = 2 and KB = 1
      (18) Same as (8), but with KA = 2 and KB = 2
      (19) Same as (8), but with KA = 3 and KB = 1
      (20) Same as (8), but with KA = 3 and KB = 2
      (21) Same as (8), but with KA = 3 and KB = 3
  NSIZES  INTEGER
          The number of sizes of matrices to use.  If it is zero,
          ZDRVSG2STG 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, ZDRVSG2STG
          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 ZDRVSG2STG 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.

  B       COMPLEX*16 array, dimension (LDB , max(NN))
          Used to hold the Hermitian positive definite matrix for
          the generailzed problem.
          On exit, B contains the last matrix actually
          used.
          Modified.

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

  D       DOUBLE PRECISION array, dimension (max(NN))
          The eigenvalues of A. On exit, the eigenvalues in D
          correspond with the matrix in A.
          Modified.

  Z       COMPLEX*16 array, dimension (LDZ, max(NN))
          The matrix of eigenvectors.
          Modified.

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

  AB      COMPLEX*16 array, dimension (LDA, max(NN))
          Workspace.
          Modified.

  BB      COMPLEX*16 array, dimension (LDB, max(NN))
          Workspace.
          Modified.

  AP      COMPLEX*16 array, dimension (max(NN)**2)
          Workspace.
          Modified.

  BP      COMPLEX*16 array, dimension (max(NN)**2)
          Workspace.
          Modified.

  WORK    COMPLEX*16 array, dimension (NWORK)
          Workspace.
          Modified.

  NWORK   INTEGER
          The number of entries in WORK.  This must be at least
          2*N + N**2  where  N = max( NN(j), 2 ).
          Not modified.

  RWORK   DOUBLE PRECISION array, dimension (LRWORK)
          Workspace.
          Modified.

  LRWORK  INTEGER
          The number of entries in RWORK.  This must be at least
          max( 7*N, 1 + 4*N + 2*N*lg(N) + 3*N**2 ) where
          N = max( NN(j) ) and lg( N ) = smallest integer k such
          that 2**k >= N .
          Not modified.

  IWORK   INTEGER array, dimension (LIWORK))
          Workspace.
          Modified.

  LIWORK  INTEGER
          The number of entries in IWORK.  This must be at least
          2 + 5*max( NN(j) ).
          Not modified.

  RESULT  DOUBLE PRECISION array, dimension (70)
          The values computed by the 70 tests described above.
          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: LDZ < 1 or LDZ < NMAX.
          -21: NWORK too small.
          -23: LRWORK too small.
          -25: LIWORK too small.
          If  ZLATMR, CLATMS, ZHEGV, ZHPGV, ZHBGV, CHEGVD, CHPGVD,
              ZHPGVD, ZHEGVX, CHPGVX, ZHBGVX 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 that have been run
                       on this matrix.
       NTESTT          The total number of tests for this call.
       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 372 of file zdrvsg2stg.f.

376*
377 IMPLICIT NONE
378*
379* -- LAPACK test routine --
380* -- LAPACK is a software package provided by Univ. of Tennessee, --
381* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
382*
383* .. Scalar Arguments ..
384 INTEGER INFO, LDA, LDB, LDZ, LIWORK, LRWORK, NOUNIT,
385 $ NSIZES, NTYPES, NWORK
386 DOUBLE PRECISION THRESH
387* ..
388* .. Array Arguments ..
389 LOGICAL DOTYPE( * )
390 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
391 DOUBLE PRECISION D( * ), D2( * ), RESULT( * ), RWORK( * )
392 COMPLEX*16 A( LDA, * ), AB( LDA, * ), AP( * ),
393 $ B( LDB, * ), BB( LDB, * ), BP( * ), WORK( * ),
394 $ Z( LDZ, * )
395* ..
396*
397* =====================================================================
398*
399* .. Parameters ..
400 DOUBLE PRECISION ZERO, ONE, TEN
401 parameter( zero = 0.0d+0, one = 1.0d+0, ten = 10.0d+0 )
402 COMPLEX*16 CZERO, CONE
403 parameter( czero = ( 0.0d+0, 0.0d+0 ),
404 $ cone = ( 1.0d+0, 0.0d+0 ) )
405 INTEGER MAXTYP
406 parameter( maxtyp = 21 )
407* ..
408* .. Local Scalars ..
409 LOGICAL BADNN
410 CHARACTER UPLO
411 INTEGER I, IBTYPE, IBUPLO, IINFO, IJ, IL, IMODE, ITEMP,
412 $ ITYPE, IU, J, JCOL, JSIZE, JTYPE, KA, KA9, KB,
413 $ KB9, M, MTYPES, N, NERRS, NMATS, NMAX, NTEST,
414 $ NTESTT
415 DOUBLE PRECISION ABSTOL, ANINV, ANORM, COND, OVFL, RTOVFL,
416 $ RTUNFL, ULP, ULPINV, UNFL, VL, VU, TEMP1, TEMP2
417* ..
418* .. Local Arrays ..
419 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), ISEED2( 4 ),
420 $ KMAGN( MAXTYP ), KMODE( MAXTYP ),
421 $ KTYPE( MAXTYP )
422* ..
423* .. External Functions ..
424 LOGICAL LSAME
425 DOUBLE PRECISION DLAMCH, DLARND
426 EXTERNAL lsame, dlamch, dlarnd
427* ..
428* .. External Subroutines ..
429 EXTERNAL dlabad, dlafts, dlasum, xerbla, zhbgv, zhbgvd,
433* ..
434* .. Intrinsic Functions ..
435 INTRINSIC abs, dble, max, min, sqrt
436* ..
437* .. Data statements ..
438 DATA ktype / 1, 2, 5*4, 5*5, 3*8, 6*9 /
439 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
440 $ 2, 3, 6*1 /
441 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
442 $ 0, 0, 6*4 /
443* ..
444* .. Executable Statements ..
445*
446* 1) Check for errors
447*
448 ntestt = 0
449 info = 0
450*
451 badnn = .false.
452 nmax = 0
453 DO 10 j = 1, nsizes
454 nmax = max( nmax, nn( j ) )
455 IF( nn( j ).LT.0 )
456 $ badnn = .true.
457 10 CONTINUE
458*
459* Check for errors
460*
461 IF( nsizes.LT.0 ) THEN
462 info = -1
463 ELSE IF( badnn ) THEN
464 info = -2
465 ELSE IF( ntypes.LT.0 ) THEN
466 info = -3
467 ELSE IF( lda.LE.1 .OR. lda.LT.nmax ) THEN
468 info = -9
469 ELSE IF( ldz.LE.1 .OR. ldz.LT.nmax ) THEN
470 info = -16
471 ELSE IF( 2*max( nmax, 2 )**2.GT.nwork ) THEN
472 info = -21
473 ELSE IF( 2*max( nmax, 2 )**2.GT.lrwork ) THEN
474 info = -23
475 ELSE IF( 2*max( nmax, 2 )**2.GT.liwork ) THEN
476 info = -25
477 END IF
478*
479 IF( info.NE.0 ) THEN
480 CALL xerbla( 'ZDRVSG2STG', -info )
481 RETURN
482 END IF
483*
484* Quick return if possible
485*
486 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
487 $ RETURN
488*
489* More Important constants
490*
491 unfl = dlamch( 'Safe minimum' )
492 ovfl = dlamch( 'Overflow' )
493 CALL dlabad( unfl, ovfl )
494 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
495 ulpinv = one / ulp
496 rtunfl = sqrt( unfl )
497 rtovfl = sqrt( ovfl )
498*
499 DO 20 i = 1, 4
500 iseed2( i ) = iseed( i )
501 20 CONTINUE
502*
503* Loop over sizes, types
504*
505 nerrs = 0
506 nmats = 0
507*
508 DO 650 jsize = 1, nsizes
509 n = nn( jsize )
510 aninv = one / dble( max( 1, n ) )
511*
512 IF( nsizes.NE.1 ) THEN
513 mtypes = min( maxtyp, ntypes )
514 ELSE
515 mtypes = min( maxtyp+1, ntypes )
516 END IF
517*
518 ka9 = 0
519 kb9 = 0
520 DO 640 jtype = 1, mtypes
521 IF( .NOT.dotype( jtype ) )
522 $ GO TO 640
523 nmats = nmats + 1
524 ntest = 0
525*
526 DO 30 j = 1, 4
527 ioldsd( j ) = iseed( j )
528 30 CONTINUE
529*
530* 2) Compute "A"
531*
532* Control parameters:
533*
534* KMAGN KMODE KTYPE
535* =1 O(1) clustered 1 zero
536* =2 large clustered 2 identity
537* =3 small exponential (none)
538* =4 arithmetic diagonal, w/ eigenvalues
539* =5 random log hermitian, w/ eigenvalues
540* =6 random (none)
541* =7 random diagonal
542* =8 random hermitian
543* =9 banded, w/ eigenvalues
544*
545 IF( mtypes.GT.maxtyp )
546 $ GO TO 90
547*
548 itype = ktype( jtype )
549 imode = kmode( jtype )
550*
551* Compute norm
552*
553 GO TO ( 40, 50, 60 )kmagn( jtype )
554*
555 40 CONTINUE
556 anorm = one
557 GO TO 70
558*
559 50 CONTINUE
560 anorm = ( rtovfl*ulp )*aninv
561 GO TO 70
562*
563 60 CONTINUE
564 anorm = rtunfl*n*ulpinv
565 GO TO 70
566*
567 70 CONTINUE
568*
569 iinfo = 0
570 cond = ulpinv
571*
572* Special Matrices -- Identity & Jordan block
573*
574 IF( itype.EQ.1 ) THEN
575*
576* Zero
577*
578 ka = 0
579 kb = 0
580 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
581*
582 ELSE IF( itype.EQ.2 ) THEN
583*
584* Identity
585*
586 ka = 0
587 kb = 0
588 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
589 DO 80 jcol = 1, n
590 a( jcol, jcol ) = anorm
591 80 CONTINUE
592*
593 ELSE IF( itype.EQ.4 ) THEN
594*
595* Diagonal Matrix, [Eigen]values Specified
596*
597 ka = 0
598 kb = 0
599 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
600 $ anorm, 0, 0, 'N', a, lda, work, iinfo )
601*
602 ELSE IF( itype.EQ.5 ) THEN
603*
604* Hermitian, eigenvalues specified
605*
606 ka = max( 0, n-1 )
607 kb = ka
608 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
609 $ anorm, n, n, 'N', a, lda, work, iinfo )
610*
611 ELSE IF( itype.EQ.7 ) THEN
612*
613* Diagonal, random eigenvalues
614*
615 ka = 0
616 kb = 0
617 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
618 $ 'T', 'N', work( n+1 ), 1, one,
619 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
620 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
621*
622 ELSE IF( itype.EQ.8 ) THEN
623*
624* Hermitian, random eigenvalues
625*
626 ka = max( 0, n-1 )
627 kb = ka
628 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one, cone,
629 $ 'T', 'N', work( n+1 ), 1, one,
630 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
631 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
632*
633 ELSE IF( itype.EQ.9 ) THEN
634*
635* Hermitian banded, eigenvalues specified
636*
637* The following values are used for the half-bandwidths:
638*
639* ka = 1 kb = 1
640* ka = 2 kb = 1
641* ka = 2 kb = 2
642* ka = 3 kb = 1
643* ka = 3 kb = 2
644* ka = 3 kb = 3
645*
646 kb9 = kb9 + 1
647 IF( kb9.GT.ka9 ) THEN
648 ka9 = ka9 + 1
649 kb9 = 1
650 END IF
651 ka = max( 0, min( n-1, ka9 ) )
652 kb = max( 0, min( n-1, kb9 ) )
653 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode, cond,
654 $ anorm, ka, ka, 'N', a, lda, work, iinfo )
655*
656 ELSE
657*
658 iinfo = 1
659 END IF
660*
661 IF( iinfo.NE.0 ) THEN
662 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
663 $ ioldsd
664 info = abs( iinfo )
665 RETURN
666 END IF
667*
668 90 CONTINUE
669*
670 abstol = unfl + unfl
671 IF( n.LE.1 ) THEN
672 il = 1
673 iu = n
674 ELSE
675 il = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
676 iu = 1 + int( ( n-1 )*dlarnd( 1, iseed2 ) )
677 IF( il.GT.iu ) THEN
678 itemp = il
679 il = iu
680 iu = itemp
681 END IF
682 END IF
683*
684* 3) Call ZHEGV, ZHPGV, ZHBGV, CHEGVD, CHPGVD, CHBGVD,
685* ZHEGVX, ZHPGVX and ZHBGVX, do tests.
686*
687* loop over the three generalized problems
688* IBTYPE = 1: A*x = (lambda)*B*x
689* IBTYPE = 2: A*B*x = (lambda)*x
690* IBTYPE = 3: B*A*x = (lambda)*x
691*
692 DO 630 ibtype = 1, 3
693*
694* loop over the setting UPLO
695*
696 DO 620 ibuplo = 1, 2
697 IF( ibuplo.EQ.1 )
698 $ uplo = 'U'
699 IF( ibuplo.EQ.2 )
700 $ uplo = 'L'
701*
702* Generate random well-conditioned positive definite
703* matrix B, of bandwidth not greater than that of A.
704*
705 CALL zlatms( n, n, 'U', iseed, 'P', rwork, 5, ten,
706 $ one, kb, kb, uplo, b, ldb, work( n+1 ),
707 $ iinfo )
708*
709* Test ZHEGV
710*
711 ntest = ntest + 1
712*
713 CALL zlacpy( ' ', n, n, a, lda, z, ldz )
714 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
715*
716 CALL zhegv( ibtype, 'V', uplo, n, z, ldz, bb, ldb, d,
717 $ work, nwork, rwork, iinfo )
718 IF( iinfo.NE.0 ) THEN
719 WRITE( nounit, fmt = 9999 )'ZHEGV(V,' // uplo //
720 $ ')', iinfo, n, jtype, ioldsd
721 info = abs( iinfo )
722 IF( iinfo.LT.0 ) THEN
723 RETURN
724 ELSE
725 result( ntest ) = ulpinv
726 GO TO 100
727 END IF
728 END IF
729*
730* Do Test
731*
732 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
733 $ ldz, d, work, rwork, result( ntest ) )
734*
735* Test ZHEGV_2STAGE
736*
737 ntest = ntest + 1
738*
739 CALL zlacpy( ' ', n, n, a, lda, z, ldz )
740 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
741*
742 CALL zhegv_2stage( ibtype, 'N', uplo, n, z, ldz,
743 $ bb, ldb, d2, work, nwork, rwork,
744 $ iinfo )
745 IF( iinfo.NE.0 ) THEN
746 WRITE( nounit, fmt = 9999 )
747 $ 'ZHEGV_2STAGE(V,' // uplo //
748 $ ')', iinfo, n, jtype, ioldsd
749 info = abs( iinfo )
750 IF( iinfo.LT.0 ) THEN
751 RETURN
752 ELSE
753 result( ntest ) = ulpinv
754 GO TO 100
755 END IF
756 END IF
757*
758* Do Test
759*
760C CALL ZSGT01( IBTYPE, UPLO, N, N, A, LDA, B, LDB, Z,
761C $ LDZ, D, WORK, RWORK, RESULT( NTEST ) )
762*
763* Do Tests | D1 - D2 | / ( |D1| ulp )
764* D1 computed using the standard 1-stage reduction as reference
765* D2 computed using the 2-stage reduction
766*
767 temp1 = zero
768 temp2 = zero
769 DO 151 j = 1, n
770 temp1 = max( temp1, abs( d( j ) ),
771 $ abs( d2( j ) ) )
772 temp2 = max( temp2, abs( d( j )-d2( j ) ) )
773 151 CONTINUE
774*
775 result( ntest ) = temp2 /
776 $ max( unfl, ulp*max( temp1, temp2 ) )
777*
778* Test ZHEGVD
779*
780 ntest = ntest + 1
781*
782 CALL zlacpy( ' ', n, n, a, lda, z, ldz )
783 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
784*
785 CALL zhegvd( ibtype, 'V', uplo, n, z, ldz, bb, ldb, d,
786 $ work, nwork, rwork, lrwork, iwork,
787 $ liwork, iinfo )
788 IF( iinfo.NE.0 ) THEN
789 WRITE( nounit, fmt = 9999 )'ZHEGVD(V,' // uplo //
790 $ ')', iinfo, n, jtype, ioldsd
791 info = abs( iinfo )
792 IF( iinfo.LT.0 ) THEN
793 RETURN
794 ELSE
795 result( ntest ) = ulpinv
796 GO TO 100
797 END IF
798 END IF
799*
800* Do Test
801*
802 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
803 $ ldz, d, work, rwork, result( ntest ) )
804*
805* Test ZHEGVX
806*
807 ntest = ntest + 1
808*
809 CALL zlacpy( ' ', n, n, a, lda, ab, lda )
810 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
811*
812 CALL zhegvx( ibtype, 'V', 'A', uplo, n, ab, lda, bb,
813 $ ldb, vl, vu, il, iu, abstol, m, d, z,
814 $ ldz, work, nwork, rwork, iwork( n+1 ),
815 $ iwork, iinfo )
816 IF( iinfo.NE.0 ) THEN
817 WRITE( nounit, fmt = 9999 )'ZHEGVX(V,A' // 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 100
825 END IF
826 END IF
827*
828* Do Test
829*
830 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
831 $ ldz, d, work, rwork, result( ntest ) )
832*
833 ntest = ntest + 1
834*
835 CALL zlacpy( ' ', n, n, a, lda, ab, lda )
836 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
837*
838* since we do not know the exact eigenvalues of this
839* eigenpair, we just set VL and VU as constants.
840* It is quite possible that there are no eigenvalues
841* in this interval.
842*
843 vl = zero
844 vu = anorm
845 CALL zhegvx( ibtype, 'V', 'V', uplo, n, ab, lda, bb,
846 $ ldb, vl, vu, il, iu, abstol, m, d, z,
847 $ ldz, work, nwork, rwork, iwork( n+1 ),
848 $ iwork, iinfo )
849 IF( iinfo.NE.0 ) THEN
850 WRITE( nounit, fmt = 9999 )'ZHEGVX(V,V,' //
851 $ uplo // ')', iinfo, n, jtype, ioldsd
852 info = abs( iinfo )
853 IF( iinfo.LT.0 ) THEN
854 RETURN
855 ELSE
856 result( ntest ) = ulpinv
857 GO TO 100
858 END IF
859 END IF
860*
861* Do Test
862*
863 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
864 $ ldz, d, work, rwork, result( ntest ) )
865*
866 ntest = ntest + 1
867*
868 CALL zlacpy( ' ', n, n, a, lda, ab, lda )
869 CALL zlacpy( uplo, n, n, b, ldb, bb, ldb )
870*
871 CALL zhegvx( ibtype, 'V', 'I', uplo, n, ab, lda, bb,
872 $ ldb, vl, vu, il, iu, abstol, m, d, z,
873 $ ldz, work, nwork, rwork, iwork( n+1 ),
874 $ iwork, iinfo )
875 IF( iinfo.NE.0 ) THEN
876 WRITE( nounit, fmt = 9999 )'ZHEGVX(V,I,' //
877 $ uplo // ')', iinfo, n, jtype, ioldsd
878 info = abs( iinfo )
879 IF( iinfo.LT.0 ) THEN
880 RETURN
881 ELSE
882 result( ntest ) = ulpinv
883 GO TO 100
884 END IF
885 END IF
886*
887* Do Test
888*
889 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
890 $ ldz, d, work, rwork, result( ntest ) )
891*
892 100 CONTINUE
893*
894* Test ZHPGV
895*
896 ntest = ntest + 1
897*
898* Copy the matrices into packed storage.
899*
900 IF( lsame( uplo, 'U' ) ) THEN
901 ij = 1
902 DO 120 j = 1, n
903 DO 110 i = 1, j
904 ap( ij ) = a( i, j )
905 bp( ij ) = b( i, j )
906 ij = ij + 1
907 110 CONTINUE
908 120 CONTINUE
909 ELSE
910 ij = 1
911 DO 140 j = 1, n
912 DO 130 i = j, n
913 ap( ij ) = a( i, j )
914 bp( ij ) = b( i, j )
915 ij = ij + 1
916 130 CONTINUE
917 140 CONTINUE
918 END IF
919*
920 CALL zhpgv( ibtype, 'V', uplo, n, ap, bp, d, z, ldz,
921 $ work, rwork, iinfo )
922 IF( iinfo.NE.0 ) THEN
923 WRITE( nounit, fmt = 9999 )'ZHPGV(V,' // uplo //
924 $ ')', iinfo, n, jtype, ioldsd
925 info = abs( iinfo )
926 IF( iinfo.LT.0 ) THEN
927 RETURN
928 ELSE
929 result( ntest ) = ulpinv
930 GO TO 310
931 END IF
932 END IF
933*
934* Do Test
935*
936 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
937 $ ldz, d, work, rwork, result( ntest ) )
938*
939* Test ZHPGVD
940*
941 ntest = ntest + 1
942*
943* Copy the matrices into packed storage.
944*
945 IF( lsame( uplo, 'U' ) ) THEN
946 ij = 1
947 DO 160 j = 1, n
948 DO 150 i = 1, j
949 ap( ij ) = a( i, j )
950 bp( ij ) = b( i, j )
951 ij = ij + 1
952 150 CONTINUE
953 160 CONTINUE
954 ELSE
955 ij = 1
956 DO 180 j = 1, n
957 DO 170 i = j, n
958 ap( ij ) = a( i, j )
959 bp( ij ) = b( i, j )
960 ij = ij + 1
961 170 CONTINUE
962 180 CONTINUE
963 END IF
964*
965 CALL zhpgvd( ibtype, 'V', uplo, n, ap, bp, d, z, ldz,
966 $ work, nwork, rwork, lrwork, iwork,
967 $ liwork, iinfo )
968 IF( iinfo.NE.0 ) THEN
969 WRITE( nounit, fmt = 9999 )'ZHPGVD(V,' // uplo //
970 $ ')', iinfo, n, jtype, ioldsd
971 info = abs( iinfo )
972 IF( iinfo.LT.0 ) THEN
973 RETURN
974 ELSE
975 result( ntest ) = ulpinv
976 GO TO 310
977 END IF
978 END IF
979*
980* Do Test
981*
982 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
983 $ ldz, d, work, rwork, result( ntest ) )
984*
985* Test ZHPGVX
986*
987 ntest = ntest + 1
988*
989* Copy the matrices into packed storage.
990*
991 IF( lsame( uplo, 'U' ) ) THEN
992 ij = 1
993 DO 200 j = 1, n
994 DO 190 i = 1, j
995 ap( ij ) = a( i, j )
996 bp( ij ) = b( i, j )
997 ij = ij + 1
998 190 CONTINUE
999 200 CONTINUE
1000 ELSE
1001 ij = 1
1002 DO 220 j = 1, n
1003 DO 210 i = j, n
1004 ap( ij ) = a( i, j )
1005 bp( ij ) = b( i, j )
1006 ij = ij + 1
1007 210 CONTINUE
1008 220 CONTINUE
1009 END IF
1010*
1011 CALL zhpgvx( ibtype, 'V', 'A', uplo, n, ap, bp, vl,
1012 $ vu, il, iu, abstol, m, d, z, ldz, work,
1013 $ rwork, iwork( n+1 ), iwork, info )
1014 IF( iinfo.NE.0 ) THEN
1015 WRITE( nounit, fmt = 9999 )'ZHPGVX(V,A' // uplo //
1016 $ ')', iinfo, n, jtype, ioldsd
1017 info = abs( iinfo )
1018 IF( iinfo.LT.0 ) THEN
1019 RETURN
1020 ELSE
1021 result( ntest ) = ulpinv
1022 GO TO 310
1023 END IF
1024 END IF
1025*
1026* Do Test
1027*
1028 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
1029 $ ldz, d, work, rwork, result( ntest ) )
1030*
1031 ntest = ntest + 1
1032*
1033* Copy the matrices into packed storage.
1034*
1035 IF( lsame( uplo, 'U' ) ) THEN
1036 ij = 1
1037 DO 240 j = 1, n
1038 DO 230 i = 1, j
1039 ap( ij ) = a( i, j )
1040 bp( ij ) = b( i, j )
1041 ij = ij + 1
1042 230 CONTINUE
1043 240 CONTINUE
1044 ELSE
1045 ij = 1
1046 DO 260 j = 1, n
1047 DO 250 i = j, n
1048 ap( ij ) = a( i, j )
1049 bp( ij ) = b( i, j )
1050 ij = ij + 1
1051 250 CONTINUE
1052 260 CONTINUE
1053 END IF
1054*
1055 vl = zero
1056 vu = anorm
1057 CALL zhpgvx( ibtype, 'V', 'V', uplo, n, ap, bp, vl,
1058 $ vu, il, iu, abstol, m, d, z, ldz, work,
1059 $ rwork, iwork( n+1 ), iwork, info )
1060 IF( iinfo.NE.0 ) THEN
1061 WRITE( nounit, fmt = 9999 )'ZHPGVX(V,V' // uplo //
1062 $ ')', iinfo, n, jtype, ioldsd
1063 info = abs( iinfo )
1064 IF( iinfo.LT.0 ) THEN
1065 RETURN
1066 ELSE
1067 result( ntest ) = ulpinv
1068 GO TO 310
1069 END IF
1070 END IF
1071*
1072* Do Test
1073*
1074 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
1075 $ ldz, d, work, rwork, result( ntest ) )
1076*
1077 ntest = ntest + 1
1078*
1079* Copy the matrices into packed storage.
1080*
1081 IF( lsame( uplo, 'U' ) ) THEN
1082 ij = 1
1083 DO 280 j = 1, n
1084 DO 270 i = 1, j
1085 ap( ij ) = a( i, j )
1086 bp( ij ) = b( i, j )
1087 ij = ij + 1
1088 270 CONTINUE
1089 280 CONTINUE
1090 ELSE
1091 ij = 1
1092 DO 300 j = 1, n
1093 DO 290 i = j, n
1094 ap( ij ) = a( i, j )
1095 bp( ij ) = b( i, j )
1096 ij = ij + 1
1097 290 CONTINUE
1098 300 CONTINUE
1099 END IF
1100*
1101 CALL zhpgvx( ibtype, 'V', 'I', uplo, n, ap, bp, vl,
1102 $ vu, il, iu, abstol, m, d, z, ldz, work,
1103 $ rwork, iwork( n+1 ), iwork, info )
1104 IF( iinfo.NE.0 ) THEN
1105 WRITE( nounit, fmt = 9999 )'ZHPGVX(V,I' // uplo //
1106 $ ')', iinfo, n, jtype, ioldsd
1107 info = abs( iinfo )
1108 IF( iinfo.LT.0 ) THEN
1109 RETURN
1110 ELSE
1111 result( ntest ) = ulpinv
1112 GO TO 310
1113 END IF
1114 END IF
1115*
1116* Do Test
1117*
1118 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
1119 $ ldz, d, work, rwork, result( ntest ) )
1120*
1121 310 CONTINUE
1122*
1123 IF( ibtype.EQ.1 ) THEN
1124*
1125* TEST ZHBGV
1126*
1127 ntest = ntest + 1
1128*
1129* Copy the matrices into band storage.
1130*
1131 IF( lsame( uplo, 'U' ) ) THEN
1132 DO 340 j = 1, n
1133 DO 320 i = max( 1, j-ka ), j
1134 ab( ka+1+i-j, j ) = a( i, j )
1135 320 CONTINUE
1136 DO 330 i = max( 1, j-kb ), j
1137 bb( kb+1+i-j, j ) = b( i, j )
1138 330 CONTINUE
1139 340 CONTINUE
1140 ELSE
1141 DO 370 j = 1, n
1142 DO 350 i = j, min( n, j+ka )
1143 ab( 1+i-j, j ) = a( i, j )
1144 350 CONTINUE
1145 DO 360 i = j, min( n, j+kb )
1146 bb( 1+i-j, j ) = b( i, j )
1147 360 CONTINUE
1148 370 CONTINUE
1149 END IF
1150*
1151 CALL zhbgv( 'V', uplo, n, ka, kb, ab, lda, bb, ldb,
1152 $ d, z, ldz, work, rwork, iinfo )
1153 IF( iinfo.NE.0 ) THEN
1154 WRITE( nounit, fmt = 9999 )'ZHBGV(V,' //
1155 $ uplo // ')', iinfo, n, jtype, ioldsd
1156 info = abs( iinfo )
1157 IF( iinfo.LT.0 ) THEN
1158 RETURN
1159 ELSE
1160 result( ntest ) = ulpinv
1161 GO TO 620
1162 END IF
1163 END IF
1164*
1165* Do Test
1166*
1167 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
1168 $ ldz, d, work, rwork, result( ntest ) )
1169*
1170* TEST ZHBGVD
1171*
1172 ntest = ntest + 1
1173*
1174* Copy the matrices into band storage.
1175*
1176 IF( lsame( uplo, 'U' ) ) THEN
1177 DO 400 j = 1, n
1178 DO 380 i = max( 1, j-ka ), j
1179 ab( ka+1+i-j, j ) = a( i, j )
1180 380 CONTINUE
1181 DO 390 i = max( 1, j-kb ), j
1182 bb( kb+1+i-j, j ) = b( i, j )
1183 390 CONTINUE
1184 400 CONTINUE
1185 ELSE
1186 DO 430 j = 1, n
1187 DO 410 i = j, min( n, j+ka )
1188 ab( 1+i-j, j ) = a( i, j )
1189 410 CONTINUE
1190 DO 420 i = j, min( n, j+kb )
1191 bb( 1+i-j, j ) = b( i, j )
1192 420 CONTINUE
1193 430 CONTINUE
1194 END IF
1195*
1196 CALL zhbgvd( 'V', uplo, n, ka, kb, ab, lda, bb,
1197 $ ldb, d, z, ldz, work, nwork, rwork,
1198 $ lrwork, iwork, liwork, iinfo )
1199 IF( iinfo.NE.0 ) THEN
1200 WRITE( nounit, fmt = 9999 )'ZHBGVD(V,' //
1201 $ uplo // ')', iinfo, n, jtype, ioldsd
1202 info = abs( iinfo )
1203 IF( iinfo.LT.0 ) THEN
1204 RETURN
1205 ELSE
1206 result( ntest ) = ulpinv
1207 GO TO 620
1208 END IF
1209 END IF
1210*
1211* Do Test
1212*
1213 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
1214 $ ldz, d, work, rwork, result( ntest ) )
1215*
1216* Test ZHBGVX
1217*
1218 ntest = ntest + 1
1219*
1220* Copy the matrices into band storage.
1221*
1222 IF( lsame( uplo, 'U' ) ) THEN
1223 DO 460 j = 1, n
1224 DO 440 i = max( 1, j-ka ), j
1225 ab( ka+1+i-j, j ) = a( i, j )
1226 440 CONTINUE
1227 DO 450 i = max( 1, j-kb ), j
1228 bb( kb+1+i-j, j ) = b( i, j )
1229 450 CONTINUE
1230 460 CONTINUE
1231 ELSE
1232 DO 490 j = 1, n
1233 DO 470 i = j, min( n, j+ka )
1234 ab( 1+i-j, j ) = a( i, j )
1235 470 CONTINUE
1236 DO 480 i = j, min( n, j+kb )
1237 bb( 1+i-j, j ) = b( i, j )
1238 480 CONTINUE
1239 490 CONTINUE
1240 END IF
1241*
1242 CALL zhbgvx( 'V', 'A', uplo, n, ka, kb, ab, lda,
1243 $ bb, ldb, bp, max( 1, n ), vl, vu, il,
1244 $ iu, abstol, m, d, z, ldz, work, rwork,
1245 $ iwork( n+1 ), iwork, iinfo )
1246 IF( iinfo.NE.0 ) THEN
1247 WRITE( nounit, fmt = 9999 )'ZHBGVX(V,A' //
1248 $ uplo // ')', iinfo, n, jtype, ioldsd
1249 info = abs( iinfo )
1250 IF( iinfo.LT.0 ) THEN
1251 RETURN
1252 ELSE
1253 result( ntest ) = ulpinv
1254 GO TO 620
1255 END IF
1256 END IF
1257*
1258* Do Test
1259*
1260 CALL zsgt01( ibtype, uplo, n, n, a, lda, b, ldb, z,
1261 $ ldz, d, work, rwork, result( ntest ) )
1262*
1263 ntest = ntest + 1
1264*
1265* Copy the matrices into band storage.
1266*
1267 IF( lsame( uplo, 'U' ) ) THEN
1268 DO 520 j = 1, n
1269 DO 500 i = max( 1, j-ka ), j
1270 ab( ka+1+i-j, j ) = a( i, j )
1271 500 CONTINUE
1272 DO 510 i = max( 1, j-kb ), j
1273 bb( kb+1+i-j, j ) = b( i, j )
1274 510 CONTINUE
1275 520 CONTINUE
1276 ELSE
1277 DO 550 j = 1, n
1278 DO 530 i = j, min( n, j+ka )
1279 ab( 1+i-j, j ) = a( i, j )
1280 530 CONTINUE
1281 DO 540 i = j, min( n, j+kb )
1282 bb( 1+i-j, j ) = b( i, j )
1283 540 CONTINUE
1284 550 CONTINUE
1285 END IF
1286*
1287 vl = zero
1288 vu = anorm
1289 CALL zhbgvx( 'V', 'V', uplo, n, ka, kb, ab, lda,
1290 $ bb, ldb, bp, max( 1, n ), vl, vu, il,
1291 $ iu, abstol, m, d, z, ldz, work, rwork,
1292 $ iwork( n+1 ), iwork, iinfo )
1293 IF( iinfo.NE.0 ) THEN
1294 WRITE( nounit, fmt = 9999 )'ZHBGVX(V,V' //
1295 $ uplo // ')', iinfo, n, jtype, ioldsd
1296 info = abs( iinfo )
1297 IF( iinfo.LT.0 ) THEN
1298 RETURN
1299 ELSE
1300 result( ntest ) = ulpinv
1301 GO TO 620
1302 END IF
1303 END IF
1304*
1305* Do Test
1306*
1307 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
1308 $ ldz, d, work, rwork, result( ntest ) )
1309*
1310 ntest = ntest + 1
1311*
1312* Copy the matrices into band storage.
1313*
1314 IF( lsame( uplo, 'U' ) ) THEN
1315 DO 580 j = 1, n
1316 DO 560 i = max( 1, j-ka ), j
1317 ab( ka+1+i-j, j ) = a( i, j )
1318 560 CONTINUE
1319 DO 570 i = max( 1, j-kb ), j
1320 bb( kb+1+i-j, j ) = b( i, j )
1321 570 CONTINUE
1322 580 CONTINUE
1323 ELSE
1324 DO 610 j = 1, n
1325 DO 590 i = j, min( n, j+ka )
1326 ab( 1+i-j, j ) = a( i, j )
1327 590 CONTINUE
1328 DO 600 i = j, min( n, j+kb )
1329 bb( 1+i-j, j ) = b( i, j )
1330 600 CONTINUE
1331 610 CONTINUE
1332 END IF
1333*
1334 CALL zhbgvx( 'V', 'I', uplo, n, ka, kb, ab, lda,
1335 $ bb, ldb, bp, max( 1, n ), vl, vu, il,
1336 $ iu, abstol, m, d, z, ldz, work, rwork,
1337 $ iwork( n+1 ), iwork, iinfo )
1338 IF( iinfo.NE.0 ) THEN
1339 WRITE( nounit, fmt = 9999 )'ZHBGVX(V,I' //
1340 $ uplo // ')', iinfo, n, jtype, ioldsd
1341 info = abs( iinfo )
1342 IF( iinfo.LT.0 ) THEN
1343 RETURN
1344 ELSE
1345 result( ntest ) = ulpinv
1346 GO TO 620
1347 END IF
1348 END IF
1349*
1350* Do Test
1351*
1352 CALL zsgt01( ibtype, uplo, n, m, a, lda, b, ldb, z,
1353 $ ldz, d, work, rwork, result( ntest ) )
1354*
1355 END IF
1356*
1357 620 CONTINUE
1358 630 CONTINUE
1359*
1360* End of Loop -- Check for RESULT(j) > THRESH
1361*
1362 ntestt = ntestt + ntest
1363 CALL dlafts( 'ZSG', n, n, jtype, ntest, result, ioldsd,
1364 $ thresh, nounit, nerrs )
1365 640 CONTINUE
1366 650 CONTINUE
1367*
1368* Summary
1369*
1370 CALL dlasum( 'ZSG', nounit, nerrs, ntestt )
1371*
1372 RETURN
1373*
1374 9999 FORMAT( ' ZDRVSG2STG: ', a, ' returned INFO=', i6, '.', / 9x,
1375 $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1376*
1377* End of ZDRVSG2STG
1378*
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
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine zsgt01(ITYPE, UPLO, N, M, A, LDA, B, LDB, Z, LDZ, D, WORK, RWORK, RESULT)
ZSGT01
Definition: zsgt01.f:152
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 zhegv_2stage(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK, LWORK, RWORK, INFO)
ZHEGV_2STAGE
Definition: zhegv_2stage.f:232
subroutine zhegvd(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHEGVD
Definition: zhegvd.f:249
subroutine zhegvx(ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
ZHEGVX
Definition: zhegvx.f:307
subroutine zhegv(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK, LWORK, RWORK, INFO)
ZHEGV
Definition: zhegv.f:181
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 zhbgvd(JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHBGVD
Definition: zhbgvd.f:252
subroutine zhbgv(JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z, LDZ, WORK, RWORK, INFO)
ZHBGV
Definition: zhbgv.f:183
subroutine zhpgvx(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)
ZHPGVX
Definition: zhpgvx.f:277
subroutine zhbgvx(JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)
ZHBGVX
Definition: zhbgvx.f:300
subroutine zhpgv(ITYPE, JOBZ, UPLO, N, AP, BP, W, Z, LDZ, WORK, RWORK, INFO)
ZHPGV
Definition: zhpgv.f:165
subroutine zhpgvd(ITYPE, JOBZ, UPLO, N, AP, BP, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHPGVD
Definition: zhpgvd.f:231
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:43
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: