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

◆ ddrvsg2stg()

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

DDRVSG2STG

Purpose:
      DDRVSG2STG checks the real symmetric generalized eigenproblem
      drivers.

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

      DSYGVD, DSPGVD and DSBGVD performed the same 14 tests.

      DSYGVX, DSPGVX and DSBGVX 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,
          DDRVSG2STG 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, DDRVSG2STG
          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 DDRVSG2STG 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       DOUBLE PRECISION 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       DOUBLE PRECISION 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       DOUBLE PRECISION array, dimension (max(NN))
          The eigenvalues of A. On exit, the eigenvalues in D
          correspond with the matrix in A.
          Modified.

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

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

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

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

  WORK    DOUBLE PRECISION 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  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: LIWORK too small.
          If  DLATMR, SLATMS, DSYGV, DSPGV, DSBGV, SSYGVD, SSPGVD,
              DSBGVD, DSYGVX, DSPGVX 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 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 358 of file ddrvsg2stg.f.

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