LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dstt22()

subroutine dstt22 ( integer n,
integer m,
integer kband,
double precision, dimension( * ) ad,
double precision, dimension( * ) ae,
double precision, dimension( * ) sd,
double precision, dimension( * ) se,
double precision, dimension( ldu, * ) u,
integer ldu,
double precision, dimension( ldwork, * ) work,
integer ldwork,
double precision, dimension( 2 ) result )

DSTT22

Purpose:
!>
!> DSTT22  checks a set of M eigenvalues and eigenvectors,
!>
!>     A U = U S
!>
!> where A is symmetric tridiagonal, the columns of U are orthogonal,
!> and S is diagonal (if KBAND=0) or symmetric tridiagonal (if KBAND=1).
!> Two tests are performed:
!>
!>    RESULT(1) = | U' A U - S | / ( |A| m ulp )
!>
!>    RESULT(2) = | I - U'U | / ( m ulp )
!> 
Parameters
[in]N
!>          N is INTEGER
!>          The size of the matrix.  If it is zero, DSTT22 does nothing.
!>          It must be at least zero.
!> 
[in]M
!>          M is INTEGER
!>          The number of eigenpairs to check.  If it is zero, DSTT22
!>          does nothing.  It must be at least zero.
!> 
[in]KBAND
!>          KBAND is INTEGER
!>          The bandwidth of the matrix S.  It may only be zero or one.
!>          If zero, then S is diagonal, and SE is not referenced.  If
!>          one, then S is symmetric tri-diagonal.
!> 
[in]AD
!>          AD is DOUBLE PRECISION array, dimension (N)
!>          The diagonal of the original (unfactored) matrix A.  A is
!>          assumed to be symmetric tridiagonal.
!> 
[in]AE
!>          AE is DOUBLE PRECISION array, dimension (N)
!>          The off-diagonal of the original (unfactored) matrix A.  A
!>          is assumed to be symmetric tridiagonal.  AE(1) is ignored,
!>          AE(2) is the (1,2) and (2,1) element, etc.
!> 
[in]SD
!>          SD is DOUBLE PRECISION array, dimension (N)
!>          The diagonal of the (symmetric tri-) diagonal matrix S.
!> 
[in]SE
!>          SE is DOUBLE PRECISION array, dimension (N)
!>          The off-diagonal of the (symmetric tri-) diagonal matrix S.
!>          Not referenced if KBSND=0.  If KBAND=1, then AE(1) is
!>          ignored, SE(2) is the (1,2) and (2,1) element, etc.
!> 
[in]U
!>          U is DOUBLE PRECISION array, dimension (LDU, N)
!>          The orthogonal matrix in the decomposition.
!> 
[in]LDU
!>          LDU is INTEGER
!>          The leading dimension of U.  LDU must be at least N.
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (LDWORK, M+1)
!> 
[in]LDWORK
!>          LDWORK is INTEGER
!>          The leading dimension of WORK.  LDWORK must be at least
!>          max(1,M).
!> 
[out]RESULT
!>          RESULT is DOUBLE PRECISION array, dimension (2)
!>          The values computed by the two tests described above.  The
!>          values are currently limited to 1/ulp, to avoid overflow.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 137 of file dstt22.f.

139*
140* -- LAPACK test routine --
141* -- LAPACK is a software package provided by Univ. of Tennessee, --
142* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
143*
144* .. Scalar Arguments ..
145 INTEGER KBAND, LDU, LDWORK, M, N
146* ..
147* .. Array Arguments ..
148 DOUBLE PRECISION AD( * ), AE( * ), RESULT( 2 ), SD( * ),
149 $ SE( * ), U( LDU, * ), WORK( LDWORK, * )
150* ..
151*
152* =====================================================================
153*
154* .. Parameters ..
155 DOUBLE PRECISION ZERO, ONE
156 parameter( zero = 0.0d0, one = 1.0d0 )
157* ..
158* .. Local Scalars ..
159 INTEGER I, J, K
160 DOUBLE PRECISION ANORM, AUKJ, ULP, UNFL, WNORM
161* ..
162* .. External Functions ..
163 DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
164 EXTERNAL dlamch, dlange, dlansy
165* ..
166* .. External Subroutines ..
167 EXTERNAL dgemm
168* ..
169* .. Intrinsic Functions ..
170 INTRINSIC abs, dble, max, min
171* ..
172* .. Executable Statements ..
173*
174 result( 1 ) = zero
175 result( 2 ) = zero
176 IF( n.LE.0 .OR. m.LE.0 )
177 $ RETURN
178*
179 unfl = dlamch( 'Safe minimum' )
180 ulp = dlamch( 'Epsilon' )
181*
182* Do Test 1
183*
184* Compute the 1-norm of A.
185*
186 IF( n.GT.1 ) THEN
187 anorm = abs( ad( 1 ) ) + abs( ae( 1 ) )
188 DO 10 j = 2, n - 1
189 anorm = max( anorm, abs( ad( j ) )+abs( ae( j ) )+
190 $ abs( ae( j-1 ) ) )
191 10 CONTINUE
192 anorm = max( anorm, abs( ad( n ) )+abs( ae( n-1 ) ) )
193 ELSE
194 anorm = abs( ad( 1 ) )
195 END IF
196 anorm = max( anorm, unfl )
197*
198* Norm of U'AU - S
199*
200 DO 40 i = 1, m
201 DO 30 j = 1, m
202 work( i, j ) = zero
203 DO 20 k = 1, n
204 aukj = ad( k )*u( k, j )
205 IF( k.NE.n )
206 $ aukj = aukj + ae( k )*u( k+1, j )
207 IF( k.NE.1 )
208 $ aukj = aukj + ae( k-1 )*u( k-1, j )
209 work( i, j ) = work( i, j ) + u( k, i )*aukj
210 20 CONTINUE
211 30 CONTINUE
212 work( i, i ) = work( i, i ) - sd( i )
213 IF( kband.EQ.1 ) THEN
214 IF( i.NE.1 )
215 $ work( i, i-1 ) = work( i, i-1 ) - se( i-1 )
216 IF( i.NE.n )
217 $ work( i, i+1 ) = work( i, i+1 ) - se( i )
218 END IF
219 40 CONTINUE
220*
221 wnorm = dlansy( '1', 'L', m, work, m, work( 1, m+1 ) )
222*
223 IF( anorm.GT.wnorm ) THEN
224 result( 1 ) = ( wnorm / anorm ) / ( m*ulp )
225 ELSE
226 IF( anorm.LT.one ) THEN
227 result( 1 ) = ( min( wnorm, m*anorm ) / anorm ) / ( m*ulp )
228 ELSE
229 result( 1 ) = min( wnorm / anorm, dble( m ) ) / ( m*ulp )
230 END IF
231 END IF
232*
233* Do Test 2
234*
235* Compute U'U - I
236*
237 CALL dgemm( 'T', 'N', m, m, n, one, u, ldu, u, ldu, zero, work,
238 $ m )
239*
240 DO 50 j = 1, m
241 work( j, j ) = work( j, j ) - one
242 50 CONTINUE
243*
244 result( 2 ) = min( dble( m ), dlange( '1', m, m, work, m, work( 1,
245 $ m+1 ) ) ) / ( m*ulp )
246*
247 RETURN
248*
249* End of DSTT22
250*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:112
double precision function dlansy(norm, uplo, n, a, lda, work)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansy.f:120
Here is the call graph for this function:
Here is the caller graph for this function: