LAPACK  3.9.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.
Date
December 2016

Definition at line 135 of file zstt21.f.

135 *
136 * -- LAPACK test routine (version 3.7.0) --
137 * -- LAPACK is a software package provided by Univ. of Tennessee, --
138 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
139 * December 2016
140 *
141 * .. Scalar Arguments ..
142  INTEGER KBAND, LDU, N
143 * ..
144 * .. Array Arguments ..
145  DOUBLE PRECISION AD( * ), AE( * ), RESULT( 2 ), RWORK( * ),
146  $ SD( * ), SE( * )
147  COMPLEX*16 U( LDU, * ), WORK( * )
148 * ..
149 *
150 * =====================================================================
151 *
152 * .. Parameters ..
153  DOUBLE PRECISION ZERO, ONE
154  parameter( zero = 0.0d+0, one = 1.0d+0 )
155  COMPLEX*16 CZERO, CONE
156  parameter( czero = ( 0.0d+0, 0.0d+0 ),
157  $ cone = ( 1.0d+0, 0.0d+0 ) )
158 * ..
159 * .. Local Scalars ..
160  INTEGER J
161  DOUBLE PRECISION ANORM, TEMP1, TEMP2, ULP, UNFL, WNORM
162 * ..
163 * .. External Functions ..
164  DOUBLE PRECISION DLAMCH, ZLANGE, ZLANHE
165  EXTERNAL dlamch, zlange, zlanhe
166 * ..
167 * .. External Subroutines ..
168  EXTERNAL zgemm, zher, zher2, zlaset
169 * ..
170 * .. Intrinsic Functions ..
171  INTRINSIC abs, dble, dcmplx, max, min
172 * ..
173 * .. Executable Statements ..
174 *
175 * 1) Constants
176 *
177  result( 1 ) = zero
178  result( 2 ) = zero
179  IF( n.LE.0 )
180  $ RETURN
181 *
182  unfl = dlamch( 'Safe minimum' )
183  ulp = dlamch( 'Precision' )
184 *
185 * Do Test 1
186 *
187 * Copy A & Compute its 1-Norm:
188 *
189  CALL zlaset( 'Full', n, n, czero, czero, work, n )
190 *
191  anorm = zero
192  temp1 = zero
193 *
194  DO 10 j = 1, n - 1
195  work( ( n+1 )*( j-1 )+1 ) = ad( j )
196  work( ( n+1 )*( j-1 )+2 ) = ae( j )
197  temp2 = abs( ae( j ) )
198  anorm = max( anorm, abs( ad( j ) )+temp1+temp2 )
199  temp1 = temp2
200  10 CONTINUE
201 *
202  work( n**2 ) = ad( n )
203  anorm = max( anorm, abs( ad( n ) )+temp1, unfl )
204 *
205 * Norm of A - USU*
206 *
207  DO 20 j = 1, n
208  CALL zher( 'L', n, -sd( j ), u( 1, j ), 1, work, n )
209  20 CONTINUE
210 *
211  IF( n.GT.1 .AND. kband.EQ.1 ) THEN
212  DO 30 j = 1, n - 1
213  CALL zher2( 'L', n, -dcmplx( se( j ) ), u( 1, j ), 1,
214  $ u( 1, j+1 ), 1, work, n )
215  30 CONTINUE
216  END IF
217 *
218  wnorm = zlanhe( '1', 'L', n, work, n, rwork )
219 *
220  IF( anorm.GT.wnorm ) THEN
221  result( 1 ) = ( wnorm / anorm ) / ( n*ulp )
222  ELSE
223  IF( anorm.LT.one ) THEN
224  result( 1 ) = ( min( wnorm, n*anorm ) / anorm ) / ( n*ulp )
225  ELSE
226  result( 1 ) = min( wnorm / anorm, dble( n ) ) / ( n*ulp )
227  END IF
228  END IF
229 *
230 * Do Test 2
231 *
232 * Compute U U**H - I
233 *
234  CALL zgemm( 'N', 'C', n, n, n, cone, u, ldu, u, ldu, czero, work,
235  $ n )
236 *
237  DO 40 j = 1, n
238  work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - cone
239  40 CONTINUE
240 *
241  result( 2 ) = min( dble( n ), zlange( '1', n, n, work, n,
242  $ rwork ) ) / ( n*ulp )
243 *
244  RETURN
245 *
246 * End of ZSTT21
247 *
Here is the call graph for this function:
Here is the caller graph for this function:
zlange
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:117
zher2
subroutine zher2(UPLO, N, ALPHA, X, INCX, Y, INCY, A, LDA)
ZHER2
Definition: zher2.f:152
zgemm
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
zlanhe
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:126
zlaset
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:108
zher
subroutine zher(UPLO, N, ALPHA, X, INCX, A, LDA)
ZHER
Definition: zher.f:137
dlamch
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:70