LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvsg ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
real  THRESH,
integer  NOUNIT,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  D,
complex, dimension( ldz, * )  Z,
integer  LDZ,
complex, dimension( lda, * )  AB,
complex, dimension( ldb, * )  BB,
complex, dimension( * )  AP,
complex, dimension( * )  BP,
complex, dimension( * )  WORK,
integer  NWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
real, dimension( * )  RESULT,
integer  INFO 
)

CDRVSG

Purpose:
      CDRVSG checks the complex Hermitian generalized eigenproblem
      drivers.

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

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

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

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

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

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

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

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

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

      When CDRVSG 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) CHEGV with ITYPE = 1 and UPLO ='U':

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

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

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

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

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

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

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

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

      CHEGVD, CHPGVD and CHBGVD performed the same 14 tests.

      CHEGVX, CHPGVX and CHBGVX 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,
          CDRVSG 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, CDRVSG
          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 CDRVSG to continue the same random number
          sequence.
          Modified.

  THRESH  REAL
          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 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 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       REAL array, dimension (max(NN))
          The eigenvalues of A. On exit, the eigenvalues in D
          correspond with the matrix in A.
          Modified.

  Z       COMPLEX 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 array, dimension (LDA, max(NN))
          Workspace.
          Modified.

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

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

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

  WORK    COMPLEX 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   REAL 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  REAL 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  CLATMR, CLATMS, CHEGV, CHPGV, CHBGV, CHEGVD, CHPGVD,
              CHPGVD, CHEGVX, CHPGVX, CHBGVX 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 SLAFTS).
       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.
Date
November 2011

Definition at line 372 of file cdrvsg.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: