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

◆ zdrvsg()

subroutine zdrvsg ( 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,
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 
)

ZDRVSG

Purpose:
      ZDRVSG 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 ZDRVSG 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 )

      (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,
          ZDRVSG 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, ZDRVSG
          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 ZDRVSG 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 366 of file zdrvsg.f.

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