LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 295 of file dchksb.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: