LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsbt21.f
Go to the documentation of this file.
1*> \brief \b DSBT21
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 DSBT21( UPLO, N, KA, KS, A, LDA, D, E, U, LDU, WORK,
12* RESULT )
13*
14* .. Scalar Arguments ..
15* CHARACTER UPLO
16* INTEGER KA, KS, LDA, LDU, N
17* ..
18* .. Array Arguments ..
19* DOUBLE PRECISION A( LDA, * ), D( * ), E( * ), RESULT( 2 ),
20* $ U( LDU, * ), WORK( * )
21* ..
22*
23*
24*> \par Purpose:
25* =============
26*>
27*> \verbatim
28*>
29*> DSBT21 generally checks a decomposition of the form
30*>
31*> A = U S U**T
32*>
33*> where **T means transpose, A is symmetric banded, U is
34*> orthogonal, and S is diagonal (if KS=0) or symmetric
35*> tridiagonal (if KS=1).
36*>
37*> Specifically:
38*>
39*> RESULT(1) = | A - U S U**T | / ( |A| n ulp ) and
40*> RESULT(2) = | I - U U**T | / ( n ulp )
41*> \endverbatim
42*
43* Arguments:
44* ==========
45*
46*> \param[in] UPLO
47*> \verbatim
48*> UPLO is CHARACTER
49*> If UPLO='U', the upper triangle of A and V will be used and
50*> the (strictly) lower triangle will not be referenced.
51*> If UPLO='L', the lower triangle of A and V will be used and
52*> the (strictly) upper triangle will not be referenced.
53*> \endverbatim
54*>
55*> \param[in] N
56*> \verbatim
57*> N is INTEGER
58*> The size of the matrix. If it is zero, DSBT21 does nothing.
59*> It must be at least zero.
60*> \endverbatim
61*>
62*> \param[in] KA
63*> \verbatim
64*> KA is INTEGER
65*> The bandwidth of the matrix A. It must be at least zero. If
66*> it is larger than N-1, then max( 0, N-1 ) will be used.
67*> \endverbatim
68*>
69*> \param[in] KS
70*> \verbatim
71*> KS is INTEGER
72*> The bandwidth of the matrix S. It may only be zero or one.
73*> If zero, then S is diagonal, and E is not referenced. If
74*> one, then S is symmetric tri-diagonal.
75*> \endverbatim
76*>
77*> \param[in] A
78*> \verbatim
79*> A is DOUBLE PRECISION array, dimension (LDA, N)
80*> The original (unfactored) matrix. It is assumed to be
81*> symmetric, and only the upper (UPLO='U') or only the lower
82*> (UPLO='L') will be referenced.
83*> \endverbatim
84*>
85*> \param[in] LDA
86*> \verbatim
87*> LDA is INTEGER
88*> The leading dimension of A. It must be at least 1
89*> and at least min( KA, N-1 ).
90*> \endverbatim
91*>
92*> \param[in] D
93*> \verbatim
94*> D is DOUBLE PRECISION array, dimension (N)
95*> The diagonal of the (symmetric tri-) diagonal matrix S.
96*> \endverbatim
97*>
98*> \param[in] E
99*> \verbatim
100*> E is DOUBLE PRECISION array, dimension (N-1)
101*> The off-diagonal of the (symmetric tri-) diagonal matrix S.
102*> E(1) is the (1,2) and (2,1) element, E(2) is the (2,3) and
103*> (3,2) element, etc.
104*> Not referenced if KS=0.
105*> \endverbatim
106*>
107*> \param[in] U
108*> \verbatim
109*> U is DOUBLE PRECISION array, dimension (LDU, N)
110*> The orthogonal matrix in the decomposition, expressed as a
111*> dense matrix (i.e., not as a product of Householder
112*> transformations, Givens transformations, etc.)
113*> \endverbatim
114*>
115*> \param[in] LDU
116*> \verbatim
117*> LDU is INTEGER
118*> The leading dimension of U. LDU must be at least N and
119*> at least 1.
120*> \endverbatim
121*>
122*> \param[out] WORK
123*> \verbatim
124*> WORK is DOUBLE PRECISION array, dimension (N**2+N)
125*> \endverbatim
126*>
127*> \param[out] RESULT
128*> \verbatim
129*> RESULT is DOUBLE PRECISION array, dimension (2)
130*> The values computed by the two tests described above. The
131*> values are currently limited to 1/ulp, to avoid overflow.
132*> \endverbatim
133*
134* Authors:
135* ========
136*
137*> \author Univ. of Tennessee
138*> \author Univ. of California Berkeley
139*> \author Univ. of Colorado Denver
140*> \author NAG Ltd.
141*
142*> \ingroup double_eig
143*
144* =====================================================================
145 SUBROUTINE dsbt21( UPLO, N, KA, KS, A, LDA, D, E, U, LDU, WORK,
146 $ RESULT )
147*
148* -- LAPACK test routine --
149* -- LAPACK is a software package provided by Univ. of Tennessee, --
150* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151*
152* .. Scalar Arguments ..
153 CHARACTER UPLO
154 INTEGER KA, KS, LDA, LDU, N
155* ..
156* .. Array Arguments ..
157 DOUBLE PRECISION A( LDA, * ), D( * ), E( * ), RESULT( 2 ),
158 $ u( ldu, * ), work( * )
159* ..
160*
161* =====================================================================
162*
163* .. Parameters ..
164 DOUBLE PRECISION ZERO, ONE
165 parameter( zero = 0.0d0, one = 1.0d0 )
166* ..
167* .. Local Scalars ..
168 LOGICAL LOWER
169 CHARACTER CUPLO
170 INTEGER IKA, J, JC, JR, LW
171 DOUBLE PRECISION ANORM, ULP, UNFL, WNORM
172* ..
173* .. External Functions ..
174 LOGICAL LSAME
175 DOUBLE PRECISION DLAMCH, DLANGE, DLANSB, DLANSP
176 EXTERNAL lsame, dlamch, dlange, dlansb, dlansp
177* ..
178* .. External Subroutines ..
179 EXTERNAL dgemm, dspr, dspr2
180* ..
181* .. Intrinsic Functions ..
182 INTRINSIC dble, max, min
183* ..
184* .. Executable Statements ..
185*
186* Constants
187*
188 result( 1 ) = zero
189 result( 2 ) = zero
190 IF( n.LE.0 )
191 $ RETURN
192*
193 ika = max( 0, min( n-1, ka ) )
194 lw = ( n*( n+1 ) ) / 2
195*
196 IF( lsame( uplo, 'U' ) ) THEN
197 lower = .false.
198 cuplo = 'U'
199 ELSE
200 lower = .true.
201 cuplo = 'L'
202 END IF
203*
204 unfl = dlamch( 'Safe minimum' )
205 ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
206*
207* Some Error Checks
208*
209* Do Test 1
210*
211* Norm of A:
212*
213 anorm = max( dlansb( '1', cuplo, n, ika, a, lda, work ), unfl )
214*
215* Compute error matrix: Error = A - U S U**T
216*
217* Copy A from SB to SP storage format.
218*
219 j = 0
220 DO 50 jc = 1, n
221 IF( lower ) THEN
222 DO 10 jr = 1, min( ika+1, n+1-jc )
223 j = j + 1
224 work( j ) = a( jr, jc )
225 10 CONTINUE
226 DO 20 jr = ika + 2, n + 1 - jc
227 j = j + 1
228 work( j ) = zero
229 20 CONTINUE
230 ELSE
231 DO 30 jr = ika + 2, jc
232 j = j + 1
233 work( j ) = zero
234 30 CONTINUE
235 DO 40 jr = min( ika, jc-1 ), 0, -1
236 j = j + 1
237 work( j ) = a( ika+1-jr, jc )
238 40 CONTINUE
239 END IF
240 50 CONTINUE
241*
242 DO 60 j = 1, n
243 CALL dspr( cuplo, n, -d( j ), u( 1, j ), 1, work )
244 60 CONTINUE
245*
246 IF( n.GT.1 .AND. ks.EQ.1 ) THEN
247 DO 70 j = 1, n - 1
248 CALL dspr2( cuplo, n, -e( j ), u( 1, j ), 1, u( 1, j+1 ), 1,
249 $ work )
250 70 CONTINUE
251 END IF
252 wnorm = dlansp( '1', cuplo, n, work, work( lw+1 ) )
253*
254 IF( anorm.GT.wnorm ) THEN
255 result( 1 ) = ( wnorm / anorm ) / ( n*ulp )
256 ELSE
257 IF( anorm.LT.one ) THEN
258 result( 1 ) = ( min( wnorm, n*anorm ) / anorm ) / ( n*ulp )
259 ELSE
260 result( 1 ) = min( wnorm / anorm, dble( n ) ) / ( n*ulp )
261 END IF
262 END IF
263*
264* Do Test 2
265*
266* Compute U U**T - I
267*
268 CALL dgemm( 'N', 'C', n, n, n, one, u, ldu, u, ldu, zero, work,
269 $ n )
270*
271 DO 80 j = 1, n
272 work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - one
273 80 CONTINUE
274*
275 result( 2 ) = min( dlange( '1', n, n, work, n, work( n**2+1 ) ),
276 $ dble( n ) ) / ( n*ulp )
277*
278 RETURN
279*
280* End of DSBT21
281*
282 END
subroutine dsbt21(uplo, n, ka, ks, a, lda, d, e, u, ldu, work, result)
DSBT21
Definition dsbt21.f:147
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dspr2(uplo, n, alpha, x, incx, y, incy, ap)
DSPR2
Definition dspr2.f:142
subroutine dspr(uplo, n, alpha, x, incx, ap)
DSPR
Definition dspr.f:127