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