LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zchkhb ( integer  NSIZES,
integer, dimension( * )  NN,
integer  NWDTHS,
integer, dimension( * )  KK,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
integer  NOUNIT,
complex*16, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  SD,
double precision, dimension( * )  SE,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

ZCHKHB

Purpose:
 ZCHKHB tests the reduction of a Hermitian band matrix to tridiagonal
 from, used with the Hermitian eigenvalue problem.

 ZHBTRD factors a Hermitian band matrix A as  U S U* , where * means
 conjugate transpose, S is symmetric tridiagonal, and U is unitary.
 ZHBTRD can use either just the lower or just the upper triangle
 of A; ZCHKHB checks both cases.

 When ZCHKHB 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 hermitian banded reduction routine.  For each
 matrix, a number of tests will be performed:

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

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

 (3)     | A - V S V* | / ( |A| n ulp )  computed by ZHBTRD 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 unitary 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 unitary 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 unitary 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) Hermitian 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,
          ZCHKHB 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,
          ZCHKHB 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, ZCHKHB
          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 ZCHKHB 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 COMPLEX*16 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 ZHBTRD.
[out]SE
          SE is DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the off-diagonal of the tridiagonal matrix
          computed by ZHBTRD.
[out]U
          U is COMPLEX*16 array, dimension (LDU, max(NN))
          Used to hold the unitary matrix computed by ZHBTRD.
[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 COMPLEX*16 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]RWORK
          RWORK is DOUBLE PRECISION array
[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 300 of file zchkhb.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: