LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zget22 ( character  TRANSA,
character  TRANSE,
character  TRANSW,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( lde, * )  E,
integer  LDE,
complex*16, dimension( * )  W,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
double precision, dimension( 2 )  RESULT 
)

ZGET22

Purpose:
 ZGET22 does an eigenvector check.

 The basic test is:

    RESULT(1) = | A E  -  E W | / ( |A| |E| ulp )

 using the 1-norm.  It also tests the normalization of E:

    RESULT(2) = max | m-norm(E(j)) - 1 | / ( n ulp )
                 j

 where E(j) is the j-th eigenvector, and m-norm is the max-norm of a
 vector.  The max-norm of a complex n-vector x in this case is the
 maximum of |re(x(i)| + |im(x(i)| over i = 1, ..., n.
Parameters
[in]TRANSA
          TRANSA is CHARACTER*1
          Specifies whether or not A is transposed.
          = 'N':  No transpose
          = 'T':  Transpose
          = 'C':  Conjugate transpose
[in]TRANSE
          TRANSE is CHARACTER*1
          Specifies whether or not E is transposed.
          = 'N':  No transpose, eigenvectors are in columns of E
          = 'T':  Transpose, eigenvectors are in rows of E
          = 'C':  Conjugate transpose, eigenvectors are in rows of E
[in]TRANSW
          TRANSW is CHARACTER*1
          Specifies whether or not W is transposed.
          = 'N':  No transpose
          = 'T':  Transpose, same as TRANSW = 'N'
          = 'C':  Conjugate transpose, use -WI(j) instead of WI(j)
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The matrix whose eigenvectors are in E.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]E
          E is COMPLEX*16 array, dimension (LDE,N)
          The matrix of eigenvectors. If TRANSE = 'N', the eigenvectors
          are stored in the columns of E, if TRANSE = 'T' or 'C', the
          eigenvectors are stored in the rows of E.
[in]LDE
          LDE is INTEGER
          The leading dimension of the array E.  LDE >= max(1,N).
[in]W
          W is COMPLEX*16 array, dimension (N)
          The eigenvalues of A.
[out]WORK
          WORK is COMPLEX*16 array, dimension (N*N)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (2)
          RESULT(1) = | A E  -  E W | / ( |A| |E| ulp )
          RESULT(2) = max | m-norm(E(j)) - 1 | / ( n ulp )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 145 of file zget22.f.

145 *
146 * -- LAPACK test routine (version 3.4.0) --
147 * -- LAPACK is a software package provided by Univ. of Tennessee, --
148 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
149 * November 2011
150 *
151 * .. Scalar Arguments ..
152  CHARACTER transa, transe, transw
153  INTEGER lda, lde, n
154 * ..
155 * .. Array Arguments ..
156  DOUBLE PRECISION result( 2 ), rwork( * )
157  COMPLEX*16 a( lda, * ), e( lde, * ), w( * ), work( * )
158 * ..
159 *
160 * =====================================================================
161 *
162 * .. Parameters ..
163  DOUBLE PRECISION zero, one
164  parameter ( zero = 0.0d+0, one = 1.0d+0 )
165  COMPLEX*16 czero, cone
166  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
167  $ cone = ( 1.0d+0, 0.0d+0 ) )
168 * ..
169 * .. Local Scalars ..
170  CHARACTER norma, norme
171  INTEGER itrnse, itrnsw, j, jcol, joff, jrow, jvec
172  DOUBLE PRECISION anorm, enorm, enrmax, enrmin, errnrm, temp1,
173  $ ulp, unfl
174  COMPLEX*16 wtemp
175 * ..
176 * .. External Functions ..
177  LOGICAL lsame
178  DOUBLE PRECISION dlamch, zlange
179  EXTERNAL lsame, dlamch, zlange
180 * ..
181 * .. External Subroutines ..
182  EXTERNAL zgemm, zlaset
183 * ..
184 * .. Intrinsic Functions ..
185  INTRINSIC abs, dble, dconjg, dimag, max, min
186 * ..
187 * .. Executable Statements ..
188 *
189 * Initialize RESULT (in case N=0)
190 *
191  result( 1 ) = zero
192  result( 2 ) = zero
193  IF( n.LE.0 )
194  $ RETURN
195 *
196  unfl = dlamch( 'Safe minimum' )
197  ulp = dlamch( 'Precision' )
198 *
199  itrnse = 0
200  itrnsw = 0
201  norma = 'O'
202  norme = 'O'
203 *
204  IF( lsame( transa, 'T' ) .OR. lsame( transa, 'C' ) ) THEN
205  norma = 'I'
206  END IF
207 *
208  IF( lsame( transe, 'T' ) ) THEN
209  itrnse = 1
210  norme = 'I'
211  ELSE IF( lsame( transe, 'C' ) ) THEN
212  itrnse = 2
213  norme = 'I'
214  END IF
215 *
216  IF( lsame( transw, 'C' ) ) THEN
217  itrnsw = 1
218  END IF
219 *
220 * Normalization of E:
221 *
222  enrmin = one / ulp
223  enrmax = zero
224  IF( itrnse.EQ.0 ) THEN
225  DO 20 jvec = 1, n
226  temp1 = zero
227  DO 10 j = 1, n
228  temp1 = max( temp1, abs( dble( e( j, jvec ) ) )+
229  $ abs( dimag( e( j, jvec ) ) ) )
230  10 CONTINUE
231  enrmin = min( enrmin, temp1 )
232  enrmax = max( enrmax, temp1 )
233  20 CONTINUE
234  ELSE
235  DO 30 jvec = 1, n
236  rwork( jvec ) = zero
237  30 CONTINUE
238 *
239  DO 50 j = 1, n
240  DO 40 jvec = 1, n
241  rwork( jvec ) = max( rwork( jvec ),
242  $ abs( dble( e( jvec, j ) ) )+
243  $ abs( dimag( e( jvec, j ) ) ) )
244  40 CONTINUE
245  50 CONTINUE
246 *
247  DO 60 jvec = 1, n
248  enrmin = min( enrmin, rwork( jvec ) )
249  enrmax = max( enrmax, rwork( jvec ) )
250  60 CONTINUE
251  END IF
252 *
253 * Norm of A:
254 *
255  anorm = max( zlange( norma, n, n, a, lda, rwork ), unfl )
256 *
257 * Norm of E:
258 *
259  enorm = max( zlange( norme, n, n, e, lde, rwork ), ulp )
260 *
261 * Norm of error:
262 *
263 * Error = AE - EW
264 *
265  CALL zlaset( 'Full', n, n, czero, czero, work, n )
266 *
267  joff = 0
268  DO 100 jcol = 1, n
269  IF( itrnsw.EQ.0 ) THEN
270  wtemp = w( jcol )
271  ELSE
272  wtemp = dconjg( w( jcol ) )
273  END IF
274 *
275  IF( itrnse.EQ.0 ) THEN
276  DO 70 jrow = 1, n
277  work( joff+jrow ) = e( jrow, jcol )*wtemp
278  70 CONTINUE
279  ELSE IF( itrnse.EQ.1 ) THEN
280  DO 80 jrow = 1, n
281  work( joff+jrow ) = e( jcol, jrow )*wtemp
282  80 CONTINUE
283  ELSE
284  DO 90 jrow = 1, n
285  work( joff+jrow ) = dconjg( e( jcol, jrow ) )*wtemp
286  90 CONTINUE
287  END IF
288  joff = joff + n
289  100 CONTINUE
290 *
291  CALL zgemm( transa, transe, n, n, n, cone, a, lda, e, lde, -cone,
292  $ work, n )
293 *
294  errnrm = zlange( 'One', n, n, work, n, rwork ) / enorm
295 *
296 * Compute RESULT(1) (avoiding under/overflow)
297 *
298  IF( anorm.GT.errnrm ) THEN
299  result( 1 ) = ( errnrm / anorm ) / ulp
300  ELSE
301  IF( anorm.LT.one ) THEN
302  result( 1 ) = ( min( errnrm, anorm ) / anorm ) / ulp
303  ELSE
304  result( 1 ) = min( errnrm / anorm, one ) / ulp
305  END IF
306  END IF
307 *
308 * Compute RESULT(2) : the normalization error in E.
309 *
310  result( 2 ) = max( abs( enrmax-one ), abs( enrmin-one ) ) /
311  $ ( dble( n )*ulp )
312 *
313  RETURN
314 *
315 * End of ZGET22
316 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
logical function lde(RI, RJ, LR)
Definition: dblat2.f:2945
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
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
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
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: