LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchkhb.f
Go to the documentation of this file.
1*> \brief \b ZCHKHB
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE ZCHKHB( NSIZES, NN, NWDTHS, KK, NTYPES, DOTYPE, ISEED,
12* THRESH, NOUNIT, A, LDA, SD, SE, U, LDU, WORK,
13* LWORK, RWORK, RESULT, INFO )
14*
15* .. Scalar Arguments ..
16* INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES,
17* $ NWDTHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER ISEED( 4 ), KK( * ), NN( * )
23* DOUBLE PRECISION RESULT( * ), RWORK( * ), SD( * ), SE( * )
24* COMPLEX*16 A( LDA, * ), U( LDU, * ), WORK( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> ZCHKHB tests the reduction of a Hermitian band matrix to tridiagonal
34*> from, used with the Hermitian eigenvalue problem.
35*>
36*> ZHBTRD factors a Hermitian band matrix A as U S U* , where * means
37*> conjugate transpose, S is symmetric tridiagonal, and U is unitary.
38*> ZHBTRD can use either just the lower or just the upper triangle
39*> of A; ZCHKHB checks both cases.
40*>
41*> When ZCHKHB is called, a number of matrix "sizes" ("n's"), a number
42*> of bandwidths ("k's"), and a number of matrix "types" are
43*> specified. For each size ("n"), each bandwidth ("k") less than or
44*> equal to "n", and each type of matrix, one matrix will be generated
45*> and used to test the hermitian banded reduction routine. For each
46*> matrix, a number of tests will be performed:
47*>
48*> (1) | A - V S V* | / ( |A| n ulp ) computed by ZHBTRD with
49*> UPLO='U'
50*>
51*> (2) | I - UU* | / ( n ulp )
52*>
53*> (3) | A - V S V* | / ( |A| n ulp ) computed by ZHBTRD with
54*> UPLO='L'
55*>
56*> (4) | I - UU* | / ( n ulp )
57*>
58*> The "sizes" are specified by an array NN(1:NSIZES); the value of
59*> each element NN(j) specifies one size.
60*> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
61*> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
62*> Currently, the list of possible types is:
63*>
64*> (1) The zero matrix.
65*> (2) The identity matrix.
66*>
67*> (3) A diagonal matrix with evenly spaced entries
68*> 1, ..., ULP and random signs.
69*> (ULP = (first number larger than 1) - 1 )
70*> (4) A diagonal matrix with geometrically spaced entries
71*> 1, ..., ULP and random signs.
72*> (5) A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
73*> and random signs.
74*>
75*> (6) Same as (4), but multiplied by SQRT( overflow threshold )
76*> (7) Same as (4), but multiplied by SQRT( underflow threshold )
77*>
78*> (8) A matrix of the form U* D U, where U is unitary and
79*> D has evenly spaced entries 1, ..., ULP with random signs
80*> on the diagonal.
81*>
82*> (9) A matrix of the form U* D U, where U is unitary and
83*> D has geometrically spaced entries 1, ..., ULP with random
84*> signs on the diagonal.
85*>
86*> (10) A matrix of the form U* D U, where U is unitary and
87*> D has "clustered" entries 1, ULP,..., ULP with random
88*> signs on the diagonal.
89*>
90*> (11) Same as (8), but multiplied by SQRT( overflow threshold )
91*> (12) Same as (8), but multiplied by SQRT( underflow threshold )
92*>
93*> (13) Hermitian matrix with random entries chosen from (-1,1).
94*> (14) Same as (13), but multiplied by SQRT( overflow threshold )
95*> (15) Same as (13), but multiplied by SQRT( underflow threshold )
96*> \endverbatim
97*
98* Arguments:
99* ==========
100*
101*> \param[in] NSIZES
102*> \verbatim
103*> NSIZES is INTEGER
104*> The number of sizes of matrices to use. If it is zero,
105*> ZCHKHB does nothing. It must be at least zero.
106*> \endverbatim
107*>
108*> \param[in] NN
109*> \verbatim
110*> NN is INTEGER array, dimension (NSIZES)
111*> An array containing the sizes to be used for the matrices.
112*> Zero values will be skipped. The values must be at least
113*> zero.
114*> \endverbatim
115*>
116*> \param[in] NWDTHS
117*> \verbatim
118*> NWDTHS is INTEGER
119*> The number of bandwidths to use. If it is zero,
120*> ZCHKHB does nothing. It must be at least zero.
121*> \endverbatim
122*>
123*> \param[in] KK
124*> \verbatim
125*> KK is INTEGER array, dimension (NWDTHS)
126*> An array containing the bandwidths to be used for the band
127*> matrices. The values must be at least zero.
128*> \endverbatim
129*>
130*> \param[in] NTYPES
131*> \verbatim
132*> NTYPES is INTEGER
133*> The number of elements in DOTYPE. If it is zero, ZCHKHB
134*> does nothing. It must be at least zero. If it is MAXTYP+1
135*> and NSIZES is 1, then an additional type, MAXTYP+1 is
136*> defined, which is to use whatever matrix is in A. This
137*> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
138*> DOTYPE(MAXTYP+1) is .TRUE. .
139*> \endverbatim
140*>
141*> \param[in] DOTYPE
142*> \verbatim
143*> DOTYPE is LOGICAL array, dimension (NTYPES)
144*> If DOTYPE(j) is .TRUE., then for each size in NN a
145*> matrix of that size and of type j will be generated.
146*> If NTYPES is smaller than the maximum number of types
147*> defined (PARAMETER MAXTYP), then types NTYPES+1 through
148*> MAXTYP will not be generated. If NTYPES is larger
149*> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
150*> will be ignored.
151*> \endverbatim
152*>
153*> \param[in,out] ISEED
154*> \verbatim
155*> ISEED is INTEGER array, dimension (4)
156*> On entry ISEED specifies the seed of the random number
157*> generator. The array elements should be between 0 and 4095;
158*> if not they will be reduced mod 4096. Also, ISEED(4) must
159*> be odd. The random number generator uses a linear
160*> congruential sequence limited to small integers, and so
161*> should produce machine independent random numbers. The
162*> values of ISEED are changed on exit, and can be used in the
163*> next call to ZCHKHB to continue the same random number
164*> sequence.
165*> \endverbatim
166*>
167*> \param[in] THRESH
168*> \verbatim
169*> THRESH is DOUBLE PRECISION
170*> A test will count as "failed" if the "error", computed as
171*> described above, exceeds THRESH. Note that the error
172*> is scaled to be O(1), so THRESH should be a reasonably
173*> small multiple of 1, e.g., 10 or 100. In particular,
174*> it should not depend on the precision (single vs. double)
175*> or the size of the matrix. It must be at least zero.
176*> \endverbatim
177*>
178*> \param[in] NOUNIT
179*> \verbatim
180*> NOUNIT is INTEGER
181*> The FORTRAN unit number for printing out error messages
182*> (e.g., if a routine returns IINFO not equal to 0.)
183*> \endverbatim
184*>
185*> \param[in,out] A
186*> \verbatim
187*> A is COMPLEX*16 array, dimension
188*> (LDA, max(NN))
189*> Used to hold the matrix whose eigenvalues are to be
190*> computed.
191*> \endverbatim
192*>
193*> \param[in] LDA
194*> \verbatim
195*> LDA is INTEGER
196*> The leading dimension of A. It must be at least 2 (not 1!)
197*> and at least max( KK )+1.
198*> \endverbatim
199*>
200*> \param[out] SD
201*> \verbatim
202*> SD is DOUBLE PRECISION array, dimension (max(NN))
203*> Used to hold the diagonal of the tridiagonal matrix computed
204*> by ZHBTRD.
205*> \endverbatim
206*>
207*> \param[out] SE
208*> \verbatim
209*> SE is DOUBLE PRECISION array, dimension (max(NN))
210*> Used to hold the off-diagonal of the tridiagonal matrix
211*> computed by ZHBTRD.
212*> \endverbatim
213*>
214*> \param[out] U
215*> \verbatim
216*> U is COMPLEX*16 array, dimension (LDU, max(NN))
217*> Used to hold the unitary matrix computed by ZHBTRD.
218*> \endverbatim
219*>
220*> \param[in] LDU
221*> \verbatim
222*> LDU is INTEGER
223*> The leading dimension of U. It must be at least 1
224*> and at least max( NN ).
225*> \endverbatim
226*>
227*> \param[out] WORK
228*> \verbatim
229*> WORK is COMPLEX*16 array, dimension (LWORK)
230*> \endverbatim
231*>
232*> \param[in] LWORK
233*> \verbatim
234*> LWORK is INTEGER
235*> The number of entries in WORK. This must be at least
236*> max( LDA+1, max(NN)+1 )*max(NN).
237*> \endverbatim
238*>
239*> \param[out] RWORK
240*> \verbatim
241*> RWORK is DOUBLE PRECISION array
242*> \endverbatim
243*>
244*> \param[out] RESULT
245*> \verbatim
246*> RESULT is DOUBLE PRECISION array, dimension (4)
247*> The values computed by the tests described above.
248*> The values are currently limited to 1/ulp, to avoid
249*> overflow.
250*> \endverbatim
251*>
252*> \param[out] INFO
253*> \verbatim
254*> INFO is INTEGER
255*> If 0, then everything ran OK.
256*>
257*>-----------------------------------------------------------------------
258*>
259*> Some Local Variables and Parameters:
260*> ---- ----- --------- --- ----------
261*> ZERO, ONE Real 0 and 1.
262*> MAXTYP The number of types defined.
263*> NTEST The number of tests performed, or which can
264*> be performed so far, for the current matrix.
265*> NTESTT The total number of tests performed so far.
266*> NMAX Largest value in NN.
267*> NMATS The number of matrices generated so far.
268*> NERRS The number of tests which have exceeded THRESH
269*> so far.
270*> COND, IMODE Values to be passed to the matrix generators.
271*> ANORM Norm of A; passed to matrix generators.
272*>
273*> OVFL, UNFL Overflow and underflow thresholds.
274*> ULP, ULPINV Finest relative precision and its inverse.
275*> RTOVFL, RTUNFL Square roots of the previous 2 values.
276*> The following four arrays decode JTYPE:
277*> KTYPE(j) The general type (1-10) for type "j".
278*> KMODE(j) The MODE value to be passed to the matrix
279*> generator for type "j".
280*> KMAGN(j) The order of magnitude ( O(1),
281*> O(overflow^(1/2) ), O(underflow^(1/2) )
282*> \endverbatim
283*
284* Authors:
285* ========
286*
287*> \author Univ. of Tennessee
288*> \author Univ. of California Berkeley
289*> \author Univ. of Colorado Denver
290*> \author NAG Ltd.
291*
292*> \ingroup complex16_eig
293*
294* =====================================================================
295 SUBROUTINE zchkhb( NSIZES, NN, NWDTHS, KK, NTYPES, DOTYPE, ISEED,
296 $ THRESH, NOUNIT, A, LDA, SD, SE, U, LDU, WORK,
297 $ LWORK, RWORK, RESULT, INFO )
298*
299* -- LAPACK test routine --
300* -- LAPACK is a software package provided by Univ. of Tennessee, --
301* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
302*
303* .. Scalar Arguments ..
304 INTEGER INFO, LDA, LDU, LWORK, NOUNIT, NSIZES, NTYPES,
305 $ NWDTHS
306 DOUBLE PRECISION THRESH
307* ..
308* .. Array Arguments ..
309 LOGICAL DOTYPE( * )
310 INTEGER ISEED( 4 ), KK( * ), NN( * )
311 DOUBLE PRECISION RESULT( * ), RWORK( * ), SD( * ), SE( * )
312 COMPLEX*16 A( LDA, * ), U( LDU, * ), WORK( * )
313* ..
314*
315* =====================================================================
316*
317* .. Parameters ..
318 COMPLEX*16 CZERO, CONE
319 PARAMETER ( CZERO = ( 0.0d+0, 0.0d+0 ),
320 $ cone = ( 1.0d+0, 0.0d+0 ) )
321 DOUBLE PRECISION ZERO, ONE, TWO, TEN
322 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0, two = 2.0d+0,
323 $ ten = 10.0d+0 )
324 DOUBLE PRECISION HALF
325 PARAMETER ( HALF = one / two )
326 INTEGER MAXTYP
327 parameter( maxtyp = 15 )
328* ..
329* .. Local Scalars ..
330 LOGICAL BADNN, BADNNB
331 INTEGER I, IINFO, IMODE, ITYPE, J, JC, JCOL, JR, JSIZE,
332 $ jtype, jwidth, k, kmax, mtypes, n, nerrs,
333 $ nmats, nmax, ntest, ntestt
334 DOUBLE PRECISION ANINV, ANORM, COND, OVFL, RTOVFL, RTUNFL,
335 $ TEMP1, ULP, ULPINV, UNFL
336* ..
337* .. Local Arrays ..
338 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KMAGN( MAXTYP ),
339 $ KMODE( MAXTYP ), KTYPE( MAXTYP )
340* ..
341* .. External Functions ..
342 DOUBLE PRECISION DLAMCH
343 EXTERNAL DLAMCH
344* ..
345* .. External Subroutines ..
346 EXTERNAL dlasum, xerbla, zhbt21, zhbtrd, zlacpy, zlaset,
347 $ zlatmr, zlatms
348* ..
349* .. Intrinsic Functions ..
350 INTRINSIC abs, dble, dconjg, max, min, sqrt
351* ..
352* .. Data statements ..
353 DATA ktype / 1, 2, 5*4, 5*5, 3*8 /
354 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
355 $ 2, 3 /
356 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
357 $ 0, 0 /
358* ..
359* .. Executable Statements ..
360*
361* Check for errors
362*
363 ntestt = 0
364 info = 0
365*
366* Important constants
367*
368 badnn = .false.
369 nmax = 1
370 DO 10 j = 1, nsizes
371 nmax = max( nmax, nn( j ) )
372 IF( nn( j ).LT.0 )
373 $ badnn = .true.
374 10 CONTINUE
375*
376 badnnb = .false.
377 kmax = 0
378 DO 20 j = 1, nsizes
379 kmax = max( kmax, kk( j ) )
380 IF( kk( j ).LT.0 )
381 $ badnnb = .true.
382 20 CONTINUE
383 kmax = min( nmax-1, kmax )
384*
385* Check for errors
386*
387 IF( nsizes.LT.0 ) THEN
388 info = -1
389 ELSE IF( badnn ) THEN
390 info = -2
391 ELSE IF( nwdths.LT.0 ) THEN
392 info = -3
393 ELSE IF( badnnb ) THEN
394 info = -4
395 ELSE IF( ntypes.LT.0 ) THEN
396 info = -5
397 ELSE IF( lda.LT.kmax+1 ) THEN
398 info = -11
399 ELSE IF( ldu.LT.nmax ) THEN
400 info = -15
401 ELSE IF( ( max( lda, nmax )+1 )*nmax.GT.lwork ) THEN
402 info = -17
403 END IF
404*
405 IF( info.NE.0 ) THEN
406 CALL xerbla( 'ZCHKHB', -info )
407 RETURN
408 END IF
409*
410* Quick return if possible
411*
412 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 .OR. nwdths.EQ.0 )
413 $ RETURN
414*
415* More Important constants
416*
417 unfl = dlamch( 'Safe minimum' )
418 ovfl = one / unfl
419 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
420 ulpinv = one / ulp
421 rtunfl = sqrt( unfl )
422 rtovfl = sqrt( ovfl )
423*
424* Loop over sizes, types
425*
426 nerrs = 0
427 nmats = 0
428*
429 DO 190 jsize = 1, nsizes
430 n = nn( jsize )
431 aninv = one / dble( max( 1, n ) )
432*
433 DO 180 jwidth = 1, nwdths
434 k = kk( jwidth )
435 IF( k.GT.n )
436 $ GO TO 180
437 k = max( 0, min( n-1, k ) )
438*
439 IF( nsizes.NE.1 ) THEN
440 mtypes = min( maxtyp, ntypes )
441 ELSE
442 mtypes = min( maxtyp+1, ntypes )
443 END IF
444*
445 DO 170 jtype = 1, mtypes
446 IF( .NOT.dotype( jtype ) )
447 $ GO TO 170
448 nmats = nmats + 1
449 ntest = 0
450*
451 DO 30 j = 1, 4
452 ioldsd( j ) = iseed( j )
453 30 CONTINUE
454*
455* Compute "A".
456* Store as "Upper"; later, we will copy to other format.
457*
458* Control parameters:
459*
460* KMAGN KMODE KTYPE
461* =1 O(1) clustered 1 zero
462* =2 large clustered 2 identity
463* =3 small exponential (none)
464* =4 arithmetic diagonal, (w/ eigenvalues)
465* =5 random log hermitian, w/ eigenvalues
466* =6 random (none)
467* =7 random diagonal
468* =8 random hermitian
469* =9 positive definite
470* =10 diagonally dominant tridiagonal
471*
472 IF( mtypes.GT.maxtyp )
473 $ GO TO 100
474*
475 itype = ktype( jtype )
476 imode = kmode( jtype )
477*
478* Compute norm
479*
480 GO TO ( 40, 50, 60 )kmagn( jtype )
481*
482 40 CONTINUE
483 anorm = one
484 GO TO 70
485*
486 50 CONTINUE
487 anorm = ( rtovfl*ulp )*aninv
488 GO TO 70
489*
490 60 CONTINUE
491 anorm = rtunfl*n*ulpinv
492 GO TO 70
493*
494 70 CONTINUE
495*
496 CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
497 iinfo = 0
498 IF( jtype.LE.15 ) THEN
499 cond = ulpinv
500 ELSE
501 cond = ulpinv*aninv / ten
502 END IF
503*
504* Special Matrices -- Identity & Jordan block
505*
506* Zero
507*
508 IF( itype.EQ.1 ) THEN
509 iinfo = 0
510*
511 ELSE IF( itype.EQ.2 ) THEN
512*
513* Identity
514*
515 DO 80 jcol = 1, n
516 a( k+1, jcol ) = anorm
517 80 CONTINUE
518*
519 ELSE IF( itype.EQ.4 ) THEN
520*
521* Diagonal Matrix, [Eigen]values Specified
522*
523 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode,
524 $ cond, anorm, 0, 0, 'Q', a( k+1, 1 ), lda,
525 $ work, iinfo )
526*
527 ELSE IF( itype.EQ.5 ) THEN
528*
529* Hermitian, eigenvalues specified
530*
531 CALL zlatms( n, n, 'S', iseed, 'H', rwork, imode,
532 $ cond, anorm, k, k, 'Q', a, lda, work,
533 $ iinfo )
534*
535 ELSE IF( itype.EQ.7 ) THEN
536*
537* Diagonal, random eigenvalues
538*
539 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one,
540 $ cone, 'T', 'N', work( n+1 ), 1, one,
541 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
542 $ zero, anorm, 'Q', a( k+1, 1 ), lda,
543 $ idumma, iinfo )
544*
545 ELSE IF( itype.EQ.8 ) THEN
546*
547* Hermitian, random eigenvalues
548*
549 CALL zlatmr( n, n, 'S', iseed, 'H', work, 6, one,
550 $ cone, 'T', 'N', work( n+1 ), 1, one,
551 $ work( 2*n+1 ), 1, one, 'N', idumma, k, k,
552 $ zero, anorm, 'Q', a, lda, idumma, iinfo )
553*
554 ELSE IF( itype.EQ.9 ) THEN
555*
556* Positive definite, eigenvalues specified.
557*
558 CALL zlatms( n, n, 'S', iseed, 'P', rwork, imode,
559 $ cond, anorm, k, k, 'Q', a, lda,
560 $ work( n+1 ), iinfo )
561*
562 ELSE IF( itype.EQ.10 ) THEN
563*
564* Positive definite tridiagonal, eigenvalues specified.
565*
566 IF( n.GT.1 )
567 $ k = max( 1, k )
568 CALL zlatms( n, n, 'S', iseed, 'P', rwork, imode,
569 $ cond, anorm, 1, 1, 'Q', a( k, 1 ), lda,
570 $ work, iinfo )
571 DO 90 i = 2, n
572 temp1 = abs( a( k, i ) ) /
573 $ sqrt( abs( a( k+1, i-1 )*a( k+1, i ) ) )
574 IF( temp1.GT.half ) THEN
575 a( k, i ) = half*sqrt( abs( a( k+1,
576 $ i-1 )*a( k+1, i ) ) )
577 END IF
578 90 CONTINUE
579*
580 ELSE
581*
582 iinfo = 1
583 END IF
584*
585 IF( iinfo.NE.0 ) THEN
586 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n,
587 $ jtype, ioldsd
588 info = abs( iinfo )
589 RETURN
590 END IF
591*
592 100 CONTINUE
593*
594* Call ZHBTRD to compute S and U from upper triangle.
595*
596 CALL zlacpy( ' ', k+1, n, a, lda, work, lda )
597*
598 ntest = 1
599 CALL zhbtrd( 'V', 'U', n, k, work, lda, sd, se, u, ldu,
600 $ work( lda*n+1 ), iinfo )
601*
602 IF( iinfo.NE.0 ) THEN
603 WRITE( nounit, fmt = 9999 )'ZHBTRD(U)', iinfo, n,
604 $ jtype, ioldsd
605 info = abs( iinfo )
606 IF( iinfo.LT.0 ) THEN
607 RETURN
608 ELSE
609 result( 1 ) = ulpinv
610 GO TO 150
611 END IF
612 END IF
613*
614* Do tests 1 and 2
615*
616 CALL zhbt21( 'Upper', n, k, 1, a, lda, sd, se, u, ldu,
617 $ work, rwork, result( 1 ) )
618*
619* Convert A from Upper-Triangle-Only storage to
620* Lower-Triangle-Only storage.
621*
622 DO 120 jc = 1, n
623 DO 110 jr = 0, min( k, n-jc )
624 a( jr+1, jc ) = dconjg( a( k+1-jr, jc+jr ) )
625 110 CONTINUE
626 120 CONTINUE
627 DO 140 jc = n + 1 - k, n
628 DO 130 jr = min( k, n-jc ) + 1, k
629 a( jr+1, jc ) = zero
630 130 CONTINUE
631 140 CONTINUE
632*
633* Call ZHBTRD to compute S and U from lower triangle
634*
635 CALL zlacpy( ' ', k+1, n, a, lda, work, lda )
636*
637 ntest = 3
638 CALL zhbtrd( 'V', 'L', n, k, work, lda, sd, se, u, ldu,
639 $ work( lda*n+1 ), iinfo )
640*
641 IF( iinfo.NE.0 ) THEN
642 WRITE( nounit, fmt = 9999 )'ZHBTRD(L)', iinfo, n,
643 $ jtype, ioldsd
644 info = abs( iinfo )
645 IF( iinfo.LT.0 ) THEN
646 RETURN
647 ELSE
648 result( 3 ) = ulpinv
649 GO TO 150
650 END IF
651 END IF
652 ntest = 4
653*
654* Do tests 3 and 4
655*
656 CALL zhbt21( 'Lower', n, k, 1, a, lda, sd, se, u, ldu,
657 $ work, rwork, result( 3 ) )
658*
659* End of Loop -- Check for RESULT(j) > THRESH
660*
661 150 CONTINUE
662 ntestt = ntestt + ntest
663*
664* Print out tests which fail.
665*
666 DO 160 jr = 1, ntest
667 IF( result( jr ).GE.thresh ) THEN
668*
669* If this is the first test to fail,
670* print a header to the data file.
671*
672 IF( nerrs.EQ.0 ) THEN
673 WRITE( nounit, fmt = 9998 )'ZHB'
674 WRITE( nounit, fmt = 9997 )
675 WRITE( nounit, fmt = 9996 )
676 WRITE( nounit, fmt = 9995 )'Hermitian'
677 WRITE( nounit, fmt = 9994 )'unitary', '*',
678 $ 'conjugate transpose', ( '*', j = 1, 4 )
679 END IF
680 nerrs = nerrs + 1
681 WRITE( nounit, fmt = 9993 )n, k, ioldsd, jtype,
682 $ jr, result( jr )
683 END IF
684 160 CONTINUE
685*
686 170 CONTINUE
687 180 CONTINUE
688 190 CONTINUE
689*
690* Summary
691*
692 CALL dlasum( 'ZHB', nounit, nerrs, ntestt )
693 RETURN
694*
695 9999 FORMAT( ' ZCHKHB: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
696 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
697 9998 FORMAT( / 1x, a3,
698 $ ' -- Complex Hermitian Banded Tridiagonal Reduction Routines'
699 $ )
700 9997 FORMAT( ' Matrix types (see DCHK23 for details): ' )
701*
702 9996 FORMAT( / ' Special Matrices:',
703 $ / ' 1=Zero matrix. ',
704 $ ' 5=Diagonal: clustered entries.',
705 $ / ' 2=Identity matrix. ',
706 $ ' 6=Diagonal: large, evenly spaced.',
707 $ / ' 3=Diagonal: evenly spaced entries. ',
708 $ ' 7=Diagonal: small, evenly spaced.',
709 $ / ' 4=Diagonal: geometr. spaced entries.' )
710 9995 FORMAT( ' Dense ', a, ' Banded Matrices:',
711 $ / ' 8=Evenly spaced eigenvals. ',
712 $ ' 12=Small, evenly spaced eigenvals.',
713 $ / ' 9=Geometrically spaced eigenvals. ',
714 $ ' 13=Matrix with random O(1) entries.',
715 $ / ' 10=Clustered eigenvalues. ',
716 $ ' 14=Matrix with large random entries.',
717 $ / ' 11=Large, evenly spaced eigenvals. ',
718 $ ' 15=Matrix with small random entries.' )
719*
720 9994 FORMAT( / ' Tests performed: (S is Tridiag, U is ', a, ',',
721 $ / 20x, a, ' means ', a, '.', / ' UPLO=''U'':',
722 $ / ' 1= | A - U S U', a1, ' | / ( |A| n ulp ) ',
723 $ ' 2= | I - U U', a1, ' | / ( n ulp )', / ' UPLO=''L'':',
724 $ / ' 3= | A - U S U', a1, ' | / ( |A| n ulp ) ',
725 $ ' 4= | I - U U', a1, ' | / ( n ulp )' )
726 9993 FORMAT( ' N=', i5, ', K=', i4, ', seed=', 4( i4, ',' ), ' type ',
727 $ i2, ', test(', i2, ')=', g10.3 )
728*
729* End of ZCHKHB
730*
731 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlasum(type, iounit, ie, nrun)
DLASUM
Definition dlasum.f:43
subroutine zhbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
ZHBTRD
Definition zhbtrd.f:163
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
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:106
subroutine zchkhb(nsizes, nn, nwdths, kk, ntypes, dotype, iseed, thresh, nounit, a, lda, sd, se, u, ldu, work, lwork, rwork, result, info)
ZCHKHB
Definition zchkhb.f:298
subroutine zhbt21(uplo, n, ka, ks, a, lda, d, e, u, ldu, work, rwork, result)
ZHBT21
Definition zhbt21.f:152
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:490
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332