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

◆ dstt21()

subroutine dstt21 ( integer  n,
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( * )  work,
double precision, dimension( 2 )  result 
)

DSTT21

Purpose:
 DSTT21 checks a decomposition of the form

    A = U S U'

 where ' means transpose, A is symmetric tridiagonal, U is orthogonal,
 and S is diagonal (if KBAND=0) or symmetric tridiagonal (if KBAND=1).
 Two tests are performed:

    RESULT(1) = | A - U S U' | / ( |A| n ulp )

    RESULT(2) = | I - UU' | / ( n ulp )
Parameters
[in]N
          N is INTEGER
          The size of the matrix.  If it is zero, DSTT21 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-1)
          The off-diagonal of the original (unfactored) matrix A.  A
          is assumed to be symmetric tridiagonal.  AE(1) is the (1,2)
          and (2,1) element, AE(2) is the (2,3) and (3,2) 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-1)
          The off-diagonal of the (symmetric tri-) diagonal matrix S.
          Not referenced if KBSND=0.  If KBAND=1, then AE(1) is the
          (1,2) and (2,1) element, SE(2) is the (2,3) and (3,2)
          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 (N*(N+1))
[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.
          RESULT(1) is always modified.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 125 of file dstt21.f.

127*
128* -- LAPACK test routine --
129* -- LAPACK is a software package provided by Univ. of Tennessee, --
130* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
131*
132* .. Scalar Arguments ..
133 INTEGER KBAND, LDU, N
134* ..
135* .. Array Arguments ..
136 DOUBLE PRECISION AD( * ), AE( * ), RESULT( 2 ), SD( * ),
137 $ SE( * ), U( LDU, * ), WORK( * )
138* ..
139*
140* =====================================================================
141*
142* .. Parameters ..
143 DOUBLE PRECISION ZERO, ONE
144 parameter( zero = 0.0d0, one = 1.0d0 )
145* ..
146* .. Local Scalars ..
147 INTEGER J
148 DOUBLE PRECISION ANORM, TEMP1, TEMP2, ULP, UNFL, WNORM
149* ..
150* .. External Functions ..
151 DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
152 EXTERNAL dlamch, dlange, dlansy
153* ..
154* .. External Subroutines ..
155 EXTERNAL dgemm, dlaset, dsyr, dsyr2
156* ..
157* .. Intrinsic Functions ..
158 INTRINSIC abs, dble, max, min
159* ..
160* .. Executable Statements ..
161*
162* 1) Constants
163*
164 result( 1 ) = zero
165 result( 2 ) = zero
166 IF( n.LE.0 )
167 $ RETURN
168*
169 unfl = dlamch( 'Safe minimum' )
170 ulp = dlamch( 'Precision' )
171*
172* Do Test 1
173*
174* Copy A & Compute its 1-Norm:
175*
176 CALL dlaset( 'Full', n, n, zero, zero, work, n )
177*
178 anorm = zero
179 temp1 = zero
180*
181 DO 10 j = 1, n - 1
182 work( ( n+1 )*( j-1 )+1 ) = ad( j )
183 work( ( n+1 )*( j-1 )+2 ) = ae( j )
184 temp2 = abs( ae( j ) )
185 anorm = max( anorm, abs( ad( j ) )+temp1+temp2 )
186 temp1 = temp2
187 10 CONTINUE
188*
189 work( n**2 ) = ad( n )
190 anorm = max( anorm, abs( ad( n ) )+temp1, unfl )
191*
192* Norm of A - USU'
193*
194 DO 20 j = 1, n
195 CALL dsyr( 'L', n, -sd( j ), u( 1, j ), 1, work, n )
196 20 CONTINUE
197*
198 IF( n.GT.1 .AND. kband.EQ.1 ) THEN
199 DO 30 j = 1, n - 1
200 CALL dsyr2( 'L', n, -se( j ), u( 1, j ), 1, u( 1, j+1 ), 1,
201 $ work, n )
202 30 CONTINUE
203 END IF
204*
205 wnorm = dlansy( '1', 'L', n, work, n, work( n**2+1 ) )
206*
207 IF( anorm.GT.wnorm ) THEN
208 result( 1 ) = ( wnorm / anorm ) / ( n*ulp )
209 ELSE
210 IF( anorm.LT.one ) THEN
211 result( 1 ) = ( min( wnorm, n*anorm ) / anorm ) / ( n*ulp )
212 ELSE
213 result( 1 ) = min( wnorm / anorm, dble( n ) ) / ( n*ulp )
214 END IF
215 END IF
216*
217* Do Test 2
218*
219* Compute UU' - I
220*
221 CALL dgemm( 'N', 'C', n, n, n, one, u, ldu, u, ldu, zero, work,
222 $ n )
223*
224 DO 40 j = 1, n
225 work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - one
226 40 CONTINUE
227*
228 result( 2 ) = min( dble( n ), dlange( '1', n, n, work, n,
229 $ work( n**2+1 ) ) ) / ( n*ulp )
230*
231 RETURN
232*
233* End of DSTT21
234*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dsyr2(uplo, n, alpha, x, incx, y, incy, a, lda)
DSYR2
Definition dsyr2.f:147
subroutine dsyr(uplo, n, alpha, x, incx, a, lda)
DSYR
Definition dsyr.f:132
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:114
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:122
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:110
Here is the call graph for this function:
Here is the caller graph for this function: