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

◆ ddrvsg()

subroutine ddrvsg ( 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( 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 
)

DDRVSG

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

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