LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dstt22.f
Go to the documentation of this file.
1 *> \brief \b DSTT22
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 DSTT22( N, M, KBAND, AD, AE, SD, SE, U, LDU, WORK,
12 * LDWORK, RESULT )
13 *
14 * .. Scalar Arguments ..
15 * INTEGER KBAND, LDU, LDWORK, M, N
16 * ..
17 * .. Array Arguments ..
18 * DOUBLE PRECISION AD( * ), AE( * ), RESULT( 2 ), SD( * ),
19 * $ SE( * ), U( LDU, * ), WORK( LDWORK, * )
20 * ..
21 *
22 *
23 *> \par Purpose:
24 * =============
25 *>
26 *> \verbatim
27 *>
28 *> DSTT22 checks a set of M eigenvalues and eigenvectors,
29 *>
30 *> A U = U S
31 *>
32 *> where A is symmetric tridiagonal, the columns of U are orthogonal,
33 *> and S is diagonal (if KBAND=0) or symmetric tridiagonal (if KBAND=1).
34 *> Two tests are performed:
35 *>
36 *> RESULT(1) = | U' A U - S | / ( |A| m ulp )
37 *>
38 *> RESULT(2) = | I - U'U | / ( m ulp )
39 *> \endverbatim
40 *
41 * Arguments:
42 * ==========
43 *
44 *> \param[in] N
45 *> \verbatim
46 *> N is INTEGER
47 *> The size of the matrix. If it is zero, DSTT22 does nothing.
48 *> It must be at least zero.
49 *> \endverbatim
50 *>
51 *> \param[in] M
52 *> \verbatim
53 *> M is INTEGER
54 *> The number of eigenpairs to check. If it is zero, DSTT22
55 *> does nothing. It must be at least zero.
56 *> \endverbatim
57 *>
58 *> \param[in] KBAND
59 *> \verbatim
60 *> KBAND is INTEGER
61 *> The bandwidth of the matrix S. It may only be zero or one.
62 *> If zero, then S is diagonal, and SE is not referenced. If
63 *> one, then S is symmetric tri-diagonal.
64 *> \endverbatim
65 *>
66 *> \param[in] AD
67 *> \verbatim
68 *> AD is DOUBLE PRECISION array, dimension (N)
69 *> The diagonal of the original (unfactored) matrix A. A is
70 *> assumed to be symmetric tridiagonal.
71 *> \endverbatim
72 *>
73 *> \param[in] AE
74 *> \verbatim
75 *> AE is DOUBLE PRECISION array, dimension (N)
76 *> The off-diagonal of the original (unfactored) matrix A. A
77 *> is assumed to be symmetric tridiagonal. AE(1) is ignored,
78 *> AE(2) is the (1,2) and (2,1) element, etc.
79 *> \endverbatim
80 *>
81 *> \param[in] SD
82 *> \verbatim
83 *> SD is DOUBLE PRECISION array, dimension (N)
84 *> The diagonal of the (symmetric tri-) diagonal matrix S.
85 *> \endverbatim
86 *>
87 *> \param[in] SE
88 *> \verbatim
89 *> SE is DOUBLE PRECISION array, dimension (N)
90 *> The off-diagonal of the (symmetric tri-) diagonal matrix S.
91 *> Not referenced if KBSND=0. If KBAND=1, then AE(1) is
92 *> ignored, SE(2) is the (1,2) and (2,1) element, etc.
93 *> \endverbatim
94 *>
95 *> \param[in] U
96 *> \verbatim
97 *> U is DOUBLE PRECISION array, dimension (LDU, N)
98 *> The orthogonal matrix in the decomposition.
99 *> \endverbatim
100 *>
101 *> \param[in] LDU
102 *> \verbatim
103 *> LDU is INTEGER
104 *> The leading dimension of U. LDU must be at least N.
105 *> \endverbatim
106 *>
107 *> \param[out] WORK
108 *> \verbatim
109 *> WORK is DOUBLE PRECISION array, dimension (LDWORK, M+1)
110 *> \endverbatim
111 *>
112 *> \param[in] LDWORK
113 *> \verbatim
114 *> LDWORK is INTEGER
115 *> The leading dimension of WORK. LDWORK must be at least
116 *> max(1,M).
117 *> \endverbatim
118 *>
119 *> \param[out] RESULT
120 *> \verbatim
121 *> RESULT is DOUBLE PRECISION array, dimension (2)
122 *> The values computed by the two tests described above. The
123 *> values are currently limited to 1/ulp, to avoid overflow.
124 *> \endverbatim
125 *
126 * Authors:
127 * ========
128 *
129 *> \author Univ. of Tennessee
130 *> \author Univ. of California Berkeley
131 *> \author Univ. of Colorado Denver
132 *> \author NAG Ltd.
133 *
134 *> \date November 2011
135 *
136 *> \ingroup double_eig
137 *
138 * =====================================================================
139  SUBROUTINE dstt22( N, M, KBAND, AD, AE, SD, SE, U, LDU, WORK,
140  $ ldwork, result )
141 *
142 * -- LAPACK test routine (version 3.4.0) --
143 * -- LAPACK is a software package provided by Univ. of Tennessee, --
144 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
145 * November 2011
146 *
147 * .. Scalar Arguments ..
148  INTEGER KBAND, LDU, LDWORK, M, N
149 * ..
150 * .. Array Arguments ..
151  DOUBLE PRECISION AD( * ), AE( * ), RESULT( 2 ), SD( * ),
152  $ se( * ), u( ldu, * ), work( ldwork, * )
153 * ..
154 *
155 * =====================================================================
156 *
157 * .. Parameters ..
158  DOUBLE PRECISION ZERO, ONE
159  parameter ( zero = 0.0d0, one = 1.0d0 )
160 * ..
161 * .. Local Scalars ..
162  INTEGER I, J, K
163  DOUBLE PRECISION ANORM, AUKJ, ULP, UNFL, WNORM
164 * ..
165 * .. External Functions ..
166  DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
167  EXTERNAL dlamch, dlange, dlansy
168 * ..
169 * .. External Subroutines ..
170  EXTERNAL dgemm
171 * ..
172 * .. Intrinsic Functions ..
173  INTRINSIC abs, dble, max, min
174 * ..
175 * .. Executable Statements ..
176 *
177  result( 1 ) = zero
178  result( 2 ) = zero
179  IF( n.LE.0 .OR. m.LE.0 )
180  $ RETURN
181 *
182  unfl = dlamch( 'Safe minimum' )
183  ulp = dlamch( 'Epsilon' )
184 *
185 * Do Test 1
186 *
187 * Compute the 1-norm of A.
188 *
189  IF( n.GT.1 ) THEN
190  anorm = abs( ad( 1 ) ) + abs( ae( 1 ) )
191  DO 10 j = 2, n - 1
192  anorm = max( anorm, abs( ad( j ) )+abs( ae( j ) )+
193  $ abs( ae( j-1 ) ) )
194  10 CONTINUE
195  anorm = max( anorm, abs( ad( n ) )+abs( ae( n-1 ) ) )
196  ELSE
197  anorm = abs( ad( 1 ) )
198  END IF
199  anorm = max( anorm, unfl )
200 *
201 * Norm of U'AU - S
202 *
203  DO 40 i = 1, m
204  DO 30 j = 1, m
205  work( i, j ) = zero
206  DO 20 k = 1, n
207  aukj = ad( k )*u( k, j )
208  IF( k.NE.n )
209  $ aukj = aukj + ae( k )*u( k+1, j )
210  IF( k.NE.1 )
211  $ aukj = aukj + ae( k-1 )*u( k-1, j )
212  work( i, j ) = work( i, j ) + u( k, i )*aukj
213  20 CONTINUE
214  30 CONTINUE
215  work( i, i ) = work( i, i ) - sd( i )
216  IF( kband.EQ.1 ) THEN
217  IF( i.NE.1 )
218  $ work( i, i-1 ) = work( i, i-1 ) - se( i-1 )
219  IF( i.NE.n )
220  $ work( i, i+1 ) = work( i, i+1 ) - se( i )
221  END IF
222  40 CONTINUE
223 *
224  wnorm = dlansy( '1', 'L', m, work, m, work( 1, m+1 ) )
225 *
226  IF( anorm.GT.wnorm ) THEN
227  result( 1 ) = ( wnorm / anorm ) / ( m*ulp )
228  ELSE
229  IF( anorm.LT.one ) THEN
230  result( 1 ) = ( min( wnorm, m*anorm ) / anorm ) / ( m*ulp )
231  ELSE
232  result( 1 ) = min( wnorm / anorm, dble( m ) ) / ( m*ulp )
233  END IF
234  END IF
235 *
236 * Do Test 2
237 *
238 * Compute U'U - I
239 *
240  CALL dgemm( 'T', 'N', m, m, n, one, u, ldu, u, ldu, zero, work,
241  $ m )
242 *
243  DO 50 j = 1, m
244  work( j, j ) = work( j, j ) - one
245  50 CONTINUE
246 *
247  result( 2 ) = min( dble( m ), dlange( '1', m, m, work, m, work( 1,
248  $ m+1 ) ) ) / ( m*ulp )
249 *
250  RETURN
251 *
252 * End of DSTT22
253 *
254  END
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dstt22(N, M, KBAND, AD, AE, SD, SE, U, LDU, WORK, LDWORK, RESULT)
DSTT22
Definition: dstt22.f:141