LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ zstt21()

subroutine zstt21 ( integer  N,
integer  KBAND,
double precision, dimension( * )  AD,
double precision, dimension( * )  AE,
double precision, dimension( * )  SD,
double precision, dimension( * )  SE,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
double precision, dimension( 2 )  RESULT 
)

ZSTT21

Purpose:
 ZSTT21  checks a decomposition of the form

    A = U S U**H

 where **H means conjugate transpose, A is real symmetric tridiagonal,
 U is unitary, and S is real and diagonal (if KBAND=0) or symmetric
 tridiagonal (if KBAND=1).  Two tests are performed:

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

    RESULT(2) = | I - U U**H | / ( n ulp )
Parameters
[in]N
          N is INTEGER
          The size of the matrix.  If it is zero, ZSTT21 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 real 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 real (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 COMPLEX*16 array, dimension (LDU, N)
          The unitary matrix in the decomposition.
[in]LDU
          LDU is INTEGER
          The leading dimension of U.  LDU must be at least N.
[out]WORK
          WORK is COMPLEX*16 array, dimension (N**2)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[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 131 of file zstt21.f.

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 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine zher(UPLO, N, ALPHA, X, INCX, A, LDA)
ZHER
Definition: zher.f:135
subroutine zher2(UPLO, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZHER2
Definition: zher2.f:150
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
double precision function zlange(NORM, M, N, A, LDA, WORK)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlange.f:115
double precision function zlanhe(NORM, UPLO, N, A, LDA, WORK)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: zlanhe.f:124
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
Here is the call graph for this function:
Here is the caller graph for this function: