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