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

◆ sdrvsg()

subroutine sdrvsg ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
real  THRESH,
integer  NOUNIT,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  D,
real, dimension( ldz, * )  Z,
integer  LDZ,
real, dimension( lda, * )  AB,
real, dimension( ldb, * )  BB,
real, dimension( * )  AP,
real, dimension( * )  BP,
real, dimension( * )  WORK,
integer  NWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
real, dimension( * )  RESULT,
integer  INFO 
)

SDRVSG

Purpose:
      SDRVSG checks the real symmetric generalized eigenproblem
      drivers.

              SSYGV computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem.

              SSYGVD computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem using a divide and conquer algorithm.

              SSYGVX computes selected eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem.

              SSPGV computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem in packed storage.

              SSPGVD computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem in packed storage using a divide and
              conquer algorithm.

              SSPGVX computes selected eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite generalized
              eigenproblem in packed storage.

              SSBGV computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite banded
              generalized eigenproblem.

              SSBGVD computes all eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite banded
              generalized eigenproblem using a divide and conquer
              algorithm.

              SSBGVX computes selected eigenvalues and, optionally,
              eigenvectors of a real symmetric-definite banded
              generalized eigenproblem.

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

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

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

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

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

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

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

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

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

      SSYGVD, SSPGVD and SSBGVD performed the same 14 tests.

      SSYGVX, SSPGVX and SSBGVX 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 orthogonal 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 orthogonal 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 orthogonal and
           D has "clustered" entries 1, ULP,..., ULP with random
           signs on the diagonal.

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

      (13) symmetric matrix with random entries chosen from (-1,1).
      (14) Same as (13), but multiplied by SQRT( overflow threshold )
      (15) Same as (13), but multiplied by SQRT( underflow threshold)

      (16) 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,
          SDRVSG 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, SDRVSG
          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 SDRVSG 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       REAL 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 and AB.  It must be at
          least 1 and at least max( NN ).
          Not modified.

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

  LDB     INTEGER
          The leading dimension of B and BB.  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       REAL array, dimension (LDZ, max(NN))
          The matrix of eigenvectors.
          Modified.

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

  AB      REAL array, dimension (LDA, max(NN))
          Workspace.
          Modified.

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

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

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

  WORK    REAL array, dimension (NWORK)
          Workspace.
          Modified.

  NWORK   INTEGER
          The number of entries in WORK.  This must be at least
          1+5*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 WORK.  This must be at least 6*N.
          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: LIWORK too small.
          If  SLATMR, SLATMS, SSYGV, SSPGV, SSBGV, SSYGVD, SSPGVD,
              SSBGVD, SSYGVX, SSPGVX or SSBGVX 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.

Definition at line 352 of file sdrvsg.f.

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