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

◆ dchksb()

subroutine dchksb ( integer  nsizes,
integer, dimension( * )  nn,
integer  nwdths,
integer, dimension( * )  kk,
integer  ntypes,
logical, dimension( * )  dotype,
integer, dimension( 4 )  iseed,
double precision  thresh,
integer  nounit,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  sd,
double precision, dimension( * )  se,
double precision, dimension( ldu, * )  u,
integer  ldu,
double precision, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  result,
integer  info 
)

DCHKSB

Purpose:
 DCHKSB tests the reduction of a symmetric band matrix to tridiagonal
 form, used with the symmetric eigenvalue problem.

 DSBTRD factors a symmetric band matrix A as  U S U' , where ' means
 transpose, S is symmetric tridiagonal, and U is orthogonal.
 DSBTRD can use either just the lower or just the upper triangle
 of A; DCHKSB checks both cases.

 When DCHKSB is called, a number of matrix "sizes" ("n's"), a number
 of bandwidths ("k's"), and a number of matrix "types" are
 specified.  For each size ("n"), each bandwidth ("k") less than or
 equal to "n", and each type of matrix, one matrix will be generated
 and used to test the symmetric banded reduction routine.  For each
 matrix, a number of tests will be performed:

 (1)     | A - V S V' | / ( |A| n ulp )  computed by DSBTRD with
                                         UPLO='U'

 (2)     | I - UU' | / ( n ulp )

 (3)     | A - V S V' | / ( |A| n ulp )  computed by DSBTRD with
                                         UPLO='L'

 (4)     | I - UU' | / ( n ulp )

 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.
 Currently, the list of possible types 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 )
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of sizes of matrices to use.  If it is zero,
          DCHKSB does nothing.  It must be at least zero.
[in]NN
          NN is 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.
[in]NWDTHS
          NWDTHS is INTEGER
          The number of bandwidths to use.  If it is zero,
          DCHKSB does nothing.  It must be at least zero.
[in]KK
          KK is INTEGER array, dimension (NWDTHS)
          An array containing the bandwidths to be used for the band
          matrices.  The values must be at least zero.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, DCHKSB
          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. .
[in]DOTYPE
          DOTYPE is 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.
[in,out]ISEED
          ISEED is 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 DCHKSB to continue the same random number
          sequence.
[in]THRESH
          THRESH is 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.
[in]NOUNIT
          NOUNIT is INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
[in,out]A
          A is DOUBLE PRECISION array, dimension
                            (LDA, max(NN))
          Used to hold the matrix whose eigenvalues are to be
          computed.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  It must be at least 2 (not 1!)
          and at least max( KK )+1.
[out]SD
          SD is DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the diagonal of the tridiagonal matrix computed
          by DSBTRD.
[out]SE
          SE is DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the off-diagonal of the tridiagonal matrix
          computed by DSBTRD.
[out]U
          U is DOUBLE PRECISION array, dimension (LDU, max(NN))
          Used to hold the orthogonal matrix computed by DSBTRD.
[in]LDU
          LDU is INTEGER
          The leading dimension of U.  It must be at least 1
          and at least max( NN ).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          max( LDA+1, max(NN)+1 )*max(NN).
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (4)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid
          overflow.
[out]INFO
          INFO is INTEGER
          If 0, then everything ran OK.

-----------------------------------------------------------------------

       Some Local Variables and Parameters:
       ---- ----- --------- --- ----------
       ZERO, ONE       Real 0 and 1.
       MAXTYP          The number of types defined.
       NTEST           The number of tests performed, or which can
                       be performed so far, for the current matrix.
       NTESTT          The total number of tests performed so far.
       NMAX            Largest value in NN.
       NMATS           The number of matrices generated so far.
       NERRS           The number of tests which have exceeded THRESH
                       so far.
       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 290 of file dchksb.f.

293*
294* -- LAPACK test routine --
295* -- LAPACK is a software package provided by Univ. of Tennessee, --
296* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
297*
298* .. Scalar Arguments ..
299 INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES,
300 $ NWDTHS
301 DOUBLE PRECISION THRESH
302* ..
303* .. Array Arguments ..
304 LOGICAL DOTYPE( * )
305 INTEGER ISEED( 4 ), KK( * ), NN( * )
306 DOUBLE PRECISION A( LDA, * ), RESULT( * ), SD( * ), SE( * ),
307 $ U( LDU, * ), WORK( * )
308* ..
309*
310* =====================================================================
311*
312* .. Parameters ..
313 DOUBLE PRECISION ZERO, ONE, TWO, TEN
314 parameter( zero = 0.0d0, one = 1.0d0, two = 2.0d0,
315 $ ten = 10.0d0 )
316 DOUBLE PRECISION HALF
317 parameter( half = one / two )
318 INTEGER MAXTYP
319 parameter( maxtyp = 15 )
320* ..
321* .. Local Scalars ..
322 LOGICAL BADNN, BADNNB
323 INTEGER I, IINFO, IMODE, ITYPE, J, JC, JCOL, JR, JSIZE,
324 $ JTYPE, JWIDTH, K, KMAX, MTYPES, N, NERRS,
325 $ NMATS, NMAX, NTEST, NTESTT
326 DOUBLE PRECISION ANINV, ANORM, COND, OVFL, RTOVFL, RTUNFL,
327 $ TEMP1, ULP, ULPINV, UNFL
328* ..
329* .. Local Arrays ..
330 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KMAGN( MAXTYP ),
331 $ KMODE( MAXTYP ), KTYPE( MAXTYP )
332* ..
333* .. External Functions ..
334 DOUBLE PRECISION DLAMCH
335 EXTERNAL dlamch
336* ..
337* .. External Subroutines ..
338 EXTERNAL dlacpy, dlaset, dlasum, dlatmr, dlatms, dsbt21,
339 $ dsbtrd, xerbla
340* ..
341* .. Intrinsic Functions ..
342 INTRINSIC abs, dble, max, min, sqrt
343* ..
344* .. Data statements ..
345 DATA ktype / 1, 2, 5*4, 5*5, 3*8 /
346 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
347 $ 2, 3 /
348 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
349 $ 0, 0 /
350* ..
351* .. Executable Statements ..
352*
353* Check for errors
354*
355 ntestt = 0
356 info = 0
357*
358* Important constants
359*
360 badnn = .false.
361 nmax = 1
362 DO 10 j = 1, nsizes
363 nmax = max( nmax, nn( j ) )
364 IF( nn( j ).LT.0 )
365 $ badnn = .true.
366 10 CONTINUE
367*
368 badnnb = .false.
369 kmax = 0
370 DO 20 j = 1, nsizes
371 kmax = max( kmax, kk( j ) )
372 IF( kk( j ).LT.0 )
373 $ badnnb = .true.
374 20 CONTINUE
375 kmax = min( nmax-1, kmax )
376*
377* Check for errors
378*
379 IF( nsizes.LT.0 ) THEN
380 info = -1
381 ELSE IF( badnn ) THEN
382 info = -2
383 ELSE IF( nwdths.LT.0 ) THEN
384 info = -3
385 ELSE IF( badnnb ) THEN
386 info = -4
387 ELSE IF( ntypes.LT.0 ) THEN
388 info = -5
389 ELSE IF( lda.LT.kmax+1 ) THEN
390 info = -11
391 ELSE IF( ldu.LT.nmax ) THEN
392 info = -15
393 ELSE IF( ( max( lda, nmax )+1 )*nmax.GT.lwork ) THEN
394 info = -17
395 END IF
396*
397 IF( info.NE.0 ) THEN
398 CALL xerbla( 'DCHKSB', -info )
399 RETURN
400 END IF
401*
402* Quick return if possible
403*
404 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 .OR. nwdths.EQ.0 )
405 $ RETURN
406*
407* More Important constants
408*
409 unfl = dlamch( 'Safe minimum' )
410 ovfl = one / unfl
411 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
412 ulpinv = one / ulp
413 rtunfl = sqrt( unfl )
414 rtovfl = sqrt( ovfl )
415*
416* Loop over sizes, types
417*
418 nerrs = 0
419 nmats = 0
420*
421 DO 190 jsize = 1, nsizes
422 n = nn( jsize )
423 aninv = one / dble( max( 1, n ) )
424*
425 DO 180 jwidth = 1, nwdths
426 k = kk( jwidth )
427 IF( k.GT.n )
428 $ GO TO 180
429 k = max( 0, min( n-1, k ) )
430*
431 IF( nsizes.NE.1 ) THEN
432 mtypes = min( maxtyp, ntypes )
433 ELSE
434 mtypes = min( maxtyp+1, ntypes )
435 END IF
436*
437 DO 170 jtype = 1, mtypes
438 IF( .NOT.dotype( jtype ) )
439 $ GO TO 170
440 nmats = nmats + 1
441 ntest = 0
442*
443 DO 30 j = 1, 4
444 ioldsd( j ) = iseed( j )
445 30 CONTINUE
446*
447* Compute "A".
448* Store as "Upper"; later, we will copy to other format.
449*
450* Control parameters:
451*
452* KMAGN KMODE KTYPE
453* =1 O(1) clustered 1 zero
454* =2 large clustered 2 identity
455* =3 small exponential (none)
456* =4 arithmetic diagonal, (w/ eigenvalues)
457* =5 random log symmetric, w/ eigenvalues
458* =6 random (none)
459* =7 random diagonal
460* =8 random symmetric
461* =9 positive definite
462* =10 diagonally dominant tridiagonal
463*
464 IF( mtypes.GT.maxtyp )
465 $ GO TO 100
466*
467 itype = ktype( jtype )
468 imode = kmode( jtype )
469*
470* Compute norm
471*
472 GO TO ( 40, 50, 60 )kmagn( jtype )
473*
474 40 CONTINUE
475 anorm = one
476 GO TO 70
477*
478 50 CONTINUE
479 anorm = ( rtovfl*ulp )*aninv
480 GO TO 70
481*
482 60 CONTINUE
483 anorm = rtunfl*n*ulpinv
484 GO TO 70
485*
486 70 CONTINUE
487*
488 CALL dlaset( 'Full', lda, n, zero, zero, a, lda )
489 iinfo = 0
490 IF( jtype.LE.15 ) THEN
491 cond = ulpinv
492 ELSE
493 cond = ulpinv*aninv / ten
494 END IF
495*
496* Special Matrices -- Identity & Jordan block
497*
498* Zero
499*
500 IF( itype.EQ.1 ) THEN
501 iinfo = 0
502*
503 ELSE IF( itype.EQ.2 ) THEN
504*
505* Identity
506*
507 DO 80 jcol = 1, n
508 a( k+1, jcol ) = anorm
509 80 CONTINUE
510*
511 ELSE IF( itype.EQ.4 ) THEN
512*
513* Diagonal Matrix, [Eigen]values Specified
514*
515 CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
516 $ anorm, 0, 0, 'Q', a( k+1, 1 ), lda,
517 $ work( n+1 ), iinfo )
518*
519 ELSE IF( itype.EQ.5 ) THEN
520*
521* Symmetric, eigenvalues specified
522*
523 CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
524 $ anorm, k, k, 'Q', a, lda, work( n+1 ),
525 $ iinfo )
526*
527 ELSE IF( itype.EQ.7 ) THEN
528*
529* Diagonal, random eigenvalues
530*
531 CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
532 $ 'T', 'N', work( n+1 ), 1, one,
533 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
534 $ zero, anorm, 'Q', a( k+1, 1 ), lda,
535 $ idumma, iinfo )
536*
537 ELSE IF( itype.EQ.8 ) THEN
538*
539* Symmetric, random eigenvalues
540*
541 CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
542 $ 'T', 'N', work( n+1 ), 1, one,
543 $ work( 2*n+1 ), 1, one, 'N', idumma, k, k,
544 $ zero, anorm, 'Q', a, lda, idumma, iinfo )
545*
546 ELSE IF( itype.EQ.9 ) THEN
547*
548* Positive definite, eigenvalues specified.
549*
550 CALL dlatms( n, n, 'S', iseed, 'P', work, imode, cond,
551 $ anorm, k, k, 'Q', a, lda, work( n+1 ),
552 $ iinfo )
553*
554 ELSE IF( itype.EQ.10 ) THEN
555*
556* Positive definite tridiagonal, eigenvalues specified.
557*
558 IF( n.GT.1 )
559 $ k = max( 1, k )
560 CALL dlatms( n, n, 'S', iseed, 'P', work, imode, cond,
561 $ anorm, 1, 1, 'Q', a( k, 1 ), lda,
562 $ work( n+1 ), iinfo )
563 DO 90 i = 2, n
564 temp1 = abs( a( k, i ) ) /
565 $ sqrt( abs( a( k+1, i-1 )*a( k+1, i ) ) )
566 IF( temp1.GT.half ) THEN
567 a( k, i ) = half*sqrt( abs( a( k+1,
568 $ i-1 )*a( k+1, i ) ) )
569 END IF
570 90 CONTINUE
571*
572 ELSE
573*
574 iinfo = 1
575 END IF
576*
577 IF( iinfo.NE.0 ) THEN
578 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n,
579 $ jtype, ioldsd
580 info = abs( iinfo )
581 RETURN
582 END IF
583*
584 100 CONTINUE
585*
586* Call DSBTRD to compute S and U from upper triangle.
587*
588 CALL dlacpy( ' ', k+1, n, a, lda, work, lda )
589*
590 ntest = 1
591 CALL dsbtrd( 'V', 'U', n, k, work, lda, sd, se, u, ldu,
592 $ work( lda*n+1 ), iinfo )
593*
594 IF( iinfo.NE.0 ) THEN
595 WRITE( nounit, fmt = 9999 )'DSBTRD(U)', iinfo, n,
596 $ jtype, ioldsd
597 info = abs( iinfo )
598 IF( iinfo.LT.0 ) THEN
599 RETURN
600 ELSE
601 result( 1 ) = ulpinv
602 GO TO 150
603 END IF
604 END IF
605*
606* Do tests 1 and 2
607*
608 CALL dsbt21( 'Upper', n, k, 1, a, lda, sd, se, u, ldu,
609 $ work, result( 1 ) )
610*
611* Convert A from Upper-Triangle-Only storage to
612* Lower-Triangle-Only storage.
613*
614 DO 120 jc = 1, n
615 DO 110 jr = 0, min( k, n-jc )
616 a( jr+1, jc ) = a( k+1-jr, jc+jr )
617 110 CONTINUE
618 120 CONTINUE
619 DO 140 jc = n + 1 - k, n
620 DO 130 jr = min( k, n-jc ) + 1, k
621 a( jr+1, jc ) = zero
622 130 CONTINUE
623 140 CONTINUE
624*
625* Call DSBTRD to compute S and U from lower triangle
626*
627 CALL dlacpy( ' ', k+1, n, a, lda, work, lda )
628*
629 ntest = 3
630 CALL dsbtrd( 'V', 'L', n, k, work, lda, sd, se, u, ldu,
631 $ work( lda*n+1 ), iinfo )
632*
633 IF( iinfo.NE.0 ) THEN
634 WRITE( nounit, fmt = 9999 )'DSBTRD(L)', iinfo, n,
635 $ jtype, ioldsd
636 info = abs( iinfo )
637 IF( iinfo.LT.0 ) THEN
638 RETURN
639 ELSE
640 result( 3 ) = ulpinv
641 GO TO 150
642 END IF
643 END IF
644 ntest = 4
645*
646* Do tests 3 and 4
647*
648 CALL dsbt21( 'Lower', n, k, 1, a, lda, sd, se, u, ldu,
649 $ work, result( 3 ) )
650*
651* End of Loop -- Check for RESULT(j) > THRESH
652*
653 150 CONTINUE
654 ntestt = ntestt + ntest
655*
656* Print out tests which fail.
657*
658 DO 160 jr = 1, ntest
659 IF( result( jr ).GE.thresh ) THEN
660*
661* If this is the first test to fail,
662* print a header to the data file.
663*
664 IF( nerrs.EQ.0 ) THEN
665 WRITE( nounit, fmt = 9998 )'DSB'
666 WRITE( nounit, fmt = 9997 )
667 WRITE( nounit, fmt = 9996 )
668 WRITE( nounit, fmt = 9995 )'Symmetric'
669 WRITE( nounit, fmt = 9994 )'orthogonal', '''',
670 $ 'transpose', ( '''', j = 1, 4 )
671 END IF
672 nerrs = nerrs + 1
673 WRITE( nounit, fmt = 9993 )n, k, ioldsd, jtype,
674 $ jr, result( jr )
675 END IF
676 160 CONTINUE
677*
678 170 CONTINUE
679 180 CONTINUE
680 190 CONTINUE
681*
682* Summary
683*
684 CALL dlasum( 'DSB', nounit, nerrs, ntestt )
685 RETURN
686*
687 9999 FORMAT( ' DCHKSB: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
688 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
689*
690 9998 FORMAT( / 1x, a3,
691 $ ' -- Real Symmetric Banded Tridiagonal Reduction Routines' )
692 9997 FORMAT( ' Matrix types (see DCHKSB for details): ' )
693*
694 9996 FORMAT( / ' Special Matrices:',
695 $ / ' 1=Zero matrix. ',
696 $ ' 5=Diagonal: clustered entries.',
697 $ / ' 2=Identity matrix. ',
698 $ ' 6=Diagonal: large, evenly spaced.',
699 $ / ' 3=Diagonal: evenly spaced entries. ',
700 $ ' 7=Diagonal: small, evenly spaced.',
701 $ / ' 4=Diagonal: geometr. spaced entries.' )
702 9995 FORMAT( ' Dense ', a, ' Banded Matrices:',
703 $ / ' 8=Evenly spaced eigenvals. ',
704 $ ' 12=Small, evenly spaced eigenvals.',
705 $ / ' 9=Geometrically spaced eigenvals. ',
706 $ ' 13=Matrix with random O(1) entries.',
707 $ / ' 10=Clustered eigenvalues. ',
708 $ ' 14=Matrix with large random entries.',
709 $ / ' 11=Large, evenly spaced eigenvals. ',
710 $ ' 15=Matrix with small random entries.' )
711*
712 9994 FORMAT( / ' Tests performed: (S is Tridiag, U is ', a, ',',
713 $ / 20x, a, ' means ', a, '.', / ' UPLO=''U'':',
714 $ / ' 1= | A - U S U', a1, ' | / ( |A| n ulp ) ',
715 $ ' 2= | I - U U', a1, ' | / ( n ulp )', / ' UPLO=''L'':',
716 $ / ' 3= | A - U S U', a1, ' | / ( |A| n ulp ) ',
717 $ ' 4= | I - U U', a1, ' | / ( n ulp )' )
718 9993 FORMAT( ' N=', i5, ', K=', i4, ', seed=', 4( i4, ',' ), ' type ',
719 $ i2, ', test(', i2, ')=', g10.3 )
720*
721* End of DCHKSB
722*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 dsbt21(uplo, n, ka, ks, a, lda, d, e, u, ldu, work, result)
DSBT21
Definition dsbt21.f:147
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:163
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
Here is the call graph for this function:
Here is the caller graph for this function: