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