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