LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dget23.f
Go to the documentation of this file.
1 *> \brief \b DGET23
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE DGET23( COMP, BALANC, JTYPE, THRESH, ISEED, NOUNIT, N,
12 * A, LDA, H, WR, WI, WR1, WI1, VL, LDVL, VR,
13 * LDVR, LRE, LDLRE, RCONDV, RCNDV1, RCDVIN,
14 * RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1, RESULT,
15 * WORK, LWORK, IWORK, INFO )
16 *
17 * .. Scalar Arguments ..
18 * LOGICAL COMP
19 * CHARACTER BALANC
20 * INTEGER INFO, JTYPE, LDA, LDLRE, LDVL, LDVR, LWORK, N,
21 * $ NOUNIT
22 * DOUBLE PRECISION THRESH
23 * ..
24 * .. Array Arguments ..
25 * INTEGER ISEED( 4 ), IWORK( * )
26 * DOUBLE PRECISION A( LDA, * ), H( LDA, * ), LRE( LDLRE, * ),
27 * $ RCDEIN( * ), RCDVIN( * ), RCNDE1( * ),
28 * $ RCNDV1( * ), RCONDE( * ), RCONDV( * ),
29 * $ RESULT( 11 ), SCALE( * ), SCALE1( * ),
30 * $ VL( LDVL, * ), VR( LDVR, * ), WI( * ),
31 * $ WI1( * ), WORK( * ), WR( * ), WR1( * )
32 * ..
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> DGET23 checks the nonsymmetric eigenvalue problem driver SGEEVX.
41 *> If COMP = .FALSE., the first 8 of the following tests will be
42 *> performed on the input matrix A, and also test 9 if LWORK is
43 *> sufficiently large.
44 *> if COMP is .TRUE. all 11 tests will be performed.
45 *>
46 *> (1) | A * VR - VR * W | / ( n |A| ulp )
47 *>
48 *> Here VR is the matrix of unit right eigenvectors.
49 *> W is a block diagonal matrix, with a 1x1 block for each
50 *> real eigenvalue and a 2x2 block for each complex conjugate
51 *> pair. If eigenvalues j and j+1 are a complex conjugate pair,
52 *> so WR(j) = WR(j+1) = wr and WI(j) = - WI(j+1) = wi, then the
53 *> 2 x 2 block corresponding to the pair will be:
54 *>
55 *> ( wr wi )
56 *> ( -wi wr )
57 *>
58 *> Such a block multiplying an n x 2 matrix ( ur ui ) on the
59 *> right will be the same as multiplying ur + i*ui by wr + i*wi.
60 *>
61 *> (2) | A**H * VL - VL * W**H | / ( n |A| ulp )
62 *>
63 *> Here VL is the matrix of unit left eigenvectors, A**H is the
64 *> conjugate transpose of A, and W is as above.
65 *>
66 *> (3) | |VR(i)| - 1 | / ulp and largest component real
67 *>
68 *> VR(i) denotes the i-th column of VR.
69 *>
70 *> (4) | |VL(i)| - 1 | / ulp and largest component real
71 *>
72 *> VL(i) denotes the i-th column of VL.
73 *>
74 *> (5) 0 if W(full) = W(partial), 1/ulp otherwise
75 *>
76 *> W(full) denotes the eigenvalues computed when VR, VL, RCONDV
77 *> and RCONDE are also computed, and W(partial) denotes the
78 *> eigenvalues computed when only some of VR, VL, RCONDV, and
79 *> RCONDE are computed.
80 *>
81 *> (6) 0 if VR(full) = VR(partial), 1/ulp otherwise
82 *>
83 *> VR(full) denotes the right eigenvectors computed when VL, RCONDV
84 *> and RCONDE are computed, and VR(partial) denotes the result
85 *> when only some of VL and RCONDV are computed.
86 *>
87 *> (7) 0 if VL(full) = VL(partial), 1/ulp otherwise
88 *>
89 *> VL(full) denotes the left eigenvectors computed when VR, RCONDV
90 *> and RCONDE are computed, and VL(partial) denotes the result
91 *> when only some of VR and RCONDV are computed.
92 *>
93 *> (8) 0 if SCALE, ILO, IHI, ABNRM (full) =
94 *> SCALE, ILO, IHI, ABNRM (partial)
95 *> 1/ulp otherwise
96 *>
97 *> SCALE, ILO, IHI and ABNRM describe how the matrix is balanced.
98 *> (full) is when VR, VL, RCONDE and RCONDV are also computed, and
99 *> (partial) is when some are not computed.
100 *>
101 *> (9) 0 if RCONDV(full) = RCONDV(partial), 1/ulp otherwise
102 *>
103 *> RCONDV(full) denotes the reciprocal condition numbers of the
104 *> right eigenvectors computed when VR, VL and RCONDE are also
105 *> computed. RCONDV(partial) denotes the reciprocal condition
106 *> numbers when only some of VR, VL and RCONDE are computed.
107 *>
108 *> (10) |RCONDV - RCDVIN| / cond(RCONDV)
109 *>
110 *> RCONDV is the reciprocal right eigenvector condition number
111 *> computed by DGEEVX and RCDVIN (the precomputed true value)
112 *> is supplied as input. cond(RCONDV) is the condition number of
113 *> RCONDV, and takes errors in computing RCONDV into account, so
114 *> that the resulting quantity should be O(ULP). cond(RCONDV) is
115 *> essentially given by norm(A)/RCONDE.
116 *>
117 *> (11) |RCONDE - RCDEIN| / cond(RCONDE)
118 *>
119 *> RCONDE is the reciprocal eigenvalue condition number
120 *> computed by DGEEVX and RCDEIN (the precomputed true value)
121 *> is supplied as input. cond(RCONDE) is the condition number
122 *> of RCONDE, and takes errors in computing RCONDE into account,
123 *> so that the resulting quantity should be O(ULP). cond(RCONDE)
124 *> is essentially given by norm(A)/RCONDV.
125 *> \endverbatim
126 *
127 * Arguments:
128 * ==========
129 *
130 *> \param[in] COMP
131 *> \verbatim
132 *> COMP is LOGICAL
133 *> COMP describes which input tests to perform:
134 *> = .FALSE. if the computed condition numbers are not to
135 *> be tested against RCDVIN and RCDEIN
136 *> = .TRUE. if they are to be compared
137 *> \endverbatim
138 *>
139 *> \param[in] BALANC
140 *> \verbatim
141 *> BALANC is CHARACTER
142 *> Describes the balancing option to be tested.
143 *> = 'N' for no permuting or diagonal scaling
144 *> = 'P' for permuting but no diagonal scaling
145 *> = 'S' for no permuting but diagonal scaling
146 *> = 'B' for permuting and diagonal scaling
147 *> \endverbatim
148 *>
149 *> \param[in] JTYPE
150 *> \verbatim
151 *> JTYPE is INTEGER
152 *> Type of input matrix. Used to label output if error occurs.
153 *> \endverbatim
154 *>
155 *> \param[in] THRESH
156 *> \verbatim
157 *> THRESH is DOUBLE PRECISION
158 *> A test will count as "failed" if the "error", computed as
159 *> described above, exceeds THRESH. Note that the error
160 *> is scaled to be O(1), so THRESH should be a reasonably
161 *> small multiple of 1, e.g., 10 or 100. In particular,
162 *> it should not depend on the precision (single vs. double)
163 *> or the size of the matrix. It must be at least zero.
164 *> \endverbatim
165 *>
166 *> \param[in] ISEED
167 *> \verbatim
168 *> ISEED is INTEGER array, dimension (4)
169 *> If COMP = .FALSE., the random number generator seed
170 *> used to produce matrix.
171 *> If COMP = .TRUE., ISEED(1) = the number of the example.
172 *> Used to label output if error occurs.
173 *> \endverbatim
174 *>
175 *> \param[in] NOUNIT
176 *> \verbatim
177 *> NOUNIT is INTEGER
178 *> The FORTRAN unit number for printing out error messages
179 *> (e.g., if a routine returns INFO not equal to 0.)
180 *> \endverbatim
181 *>
182 *> \param[in] N
183 *> \verbatim
184 *> N is INTEGER
185 *> The dimension of A. N must be at least 0.
186 *> \endverbatim
187 *>
188 *> \param[in,out] A
189 *> \verbatim
190 *> A is DOUBLE PRECISION array, dimension (LDA,N)
191 *> Used to hold the matrix whose eigenvalues are to be
192 *> computed.
193 *> \endverbatim
194 *>
195 *> \param[in] LDA
196 *> \verbatim
197 *> LDA is INTEGER
198 *> The leading dimension of A, and H. LDA must be at
199 *> least 1 and at least N.
200 *> \endverbatim
201 *>
202 *> \param[out] H
203 *> \verbatim
204 *> H is DOUBLE PRECISION array, dimension (LDA,N)
205 *> Another copy of the test matrix A, modified by DGEEVX.
206 *> \endverbatim
207 *>
208 *> \param[out] WR
209 *> \verbatim
210 *> WR is DOUBLE PRECISION array, dimension (N)
211 *> \endverbatim
212 *>
213 *> \param[out] WI
214 *> \verbatim
215 *> WI is DOUBLE PRECISION array, dimension (N)
216 *>
217 *> The real and imaginary parts of the eigenvalues of A.
218 *> On exit, WR + WI*i are the eigenvalues of the matrix in A.
219 *> \endverbatim
220 *>
221 *> \param[out] WR1
222 *> \verbatim
223 *> WR1 is DOUBLE PRECISION array, dimension (N)
224 *> \endverbatim
225 *>
226 *> \param[out] WI1
227 *> \verbatim
228 *> WI1 is DOUBLE PRECISION array, dimension (N)
229 *>
230 *> Like WR, WI, these arrays contain the eigenvalues of A,
231 *> but those computed when DGEEVX only computes a partial
232 *> eigendecomposition, i.e. not the eigenvalues and left
233 *> and right eigenvectors.
234 *> \endverbatim
235 *>
236 *> \param[out] VL
237 *> \verbatim
238 *> VL is DOUBLE PRECISION array, dimension (LDVL,N)
239 *> VL holds the computed left eigenvectors.
240 *> \endverbatim
241 *>
242 *> \param[in] LDVL
243 *> \verbatim
244 *> LDVL is INTEGER
245 *> Leading dimension of VL. Must be at least max(1,N).
246 *> \endverbatim
247 *>
248 *> \param[out] VR
249 *> \verbatim
250 *> VR is DOUBLE PRECISION array, dimension (LDVR,N)
251 *> VR holds the computed right eigenvectors.
252 *> \endverbatim
253 *>
254 *> \param[in] LDVR
255 *> \verbatim
256 *> LDVR is INTEGER
257 *> Leading dimension of VR. Must be at least max(1,N).
258 *> \endverbatim
259 *>
260 *> \param[out] LRE
261 *> \verbatim
262 *> LRE is DOUBLE PRECISION array, dimension (LDLRE,N)
263 *> LRE holds the computed right or left eigenvectors.
264 *> \endverbatim
265 *>
266 *> \param[in] LDLRE
267 *> \verbatim
268 *> LDLRE is INTEGER
269 *> Leading dimension of LRE. Must be at least max(1,N).
270 *> \endverbatim
271 *>
272 *> \param[out] RCONDV
273 *> \verbatim
274 *> RCONDV is DOUBLE PRECISION array, dimension (N)
275 *> RCONDV holds the computed reciprocal condition numbers
276 *> for eigenvectors.
277 *> \endverbatim
278 *>
279 *> \param[out] RCNDV1
280 *> \verbatim
281 *> RCNDV1 is DOUBLE PRECISION array, dimension (N)
282 *> RCNDV1 holds more computed reciprocal condition numbers
283 *> for eigenvectors.
284 *> \endverbatim
285 *>
286 *> \param[in] RCDVIN
287 *> \verbatim
288 *> RCDVIN is DOUBLE PRECISION array, dimension (N)
289 *> When COMP = .TRUE. RCDVIN holds the precomputed reciprocal
290 *> condition numbers for eigenvectors to be compared with
291 *> RCONDV.
292 *> \endverbatim
293 *>
294 *> \param[out] RCONDE
295 *> \verbatim
296 *> RCONDE is DOUBLE PRECISION array, dimension (N)
297 *> RCONDE holds the computed reciprocal condition numbers
298 *> for eigenvalues.
299 *> \endverbatim
300 *>
301 *> \param[out] RCNDE1
302 *> \verbatim
303 *> RCNDE1 is DOUBLE PRECISION array, dimension (N)
304 *> RCNDE1 holds more computed reciprocal condition numbers
305 *> for eigenvalues.
306 *> \endverbatim
307 *>
308 *> \param[in] RCDEIN
309 *> \verbatim
310 *> RCDEIN is DOUBLE PRECISION array, dimension (N)
311 *> When COMP = .TRUE. RCDEIN holds the precomputed reciprocal
312 *> condition numbers for eigenvalues to be compared with
313 *> RCONDE.
314 *> \endverbatim
315 *>
316 *> \param[out] SCALE
317 *> \verbatim
318 *> SCALE is DOUBLE PRECISION array, dimension (N)
319 *> Holds information describing balancing of matrix.
320 *> \endverbatim
321 *>
322 *> \param[out] SCALE1
323 *> \verbatim
324 *> SCALE1 is DOUBLE PRECISION array, dimension (N)
325 *> Holds information describing balancing of matrix.
326 *> \endverbatim
327 *>
328 *> \param[out] RESULT
329 *> \verbatim
330 *> RESULT is DOUBLE PRECISION array, dimension (11)
331 *> The values computed by the 11 tests described above.
332 *> The values are currently limited to 1/ulp, to avoid
333 *> overflow.
334 *> \endverbatim
335 *>
336 *> \param[out] WORK
337 *> \verbatim
338 *> WORK is DOUBLE PRECISION array, dimension (LWORK)
339 *> \endverbatim
340 *>
341 *> \param[in] LWORK
342 *> \verbatim
343 *> LWORK is INTEGER
344 *> The number of entries in WORK. This must be at least
345 *> 3*N, and 6*N+N**2 if tests 9, 10 or 11 are to be performed.
346 *> \endverbatim
347 *>
348 *> \param[out] IWORK
349 *> \verbatim
350 *> IWORK is INTEGER array, dimension (2*N)
351 *> \endverbatim
352 *>
353 *> \param[out] INFO
354 *> \verbatim
355 *> INFO is INTEGER
356 *> If 0, successful exit.
357 *> If <0, input parameter -INFO had an incorrect value.
358 *> If >0, DGEEVX returned an error code, the absolute
359 *> value of which is returned.
360 *> \endverbatim
361 *
362 * Authors:
363 * ========
364 *
365 *> \author Univ. of Tennessee
366 *> \author Univ. of California Berkeley
367 *> \author Univ. of Colorado Denver
368 *> \author NAG Ltd.
369 *
370 *> \date November 2011
371 *
372 *> \ingroup double_eig
373 *
374 * =====================================================================
375  SUBROUTINE dget23( COMP, BALANC, JTYPE, THRESH, ISEED, NOUNIT, N,
376  $ a, lda, h, wr, wi, wr1, wi1, vl, ldvl, vr,
377  $ ldvr, lre, ldlre, rcondv, rcndv1, rcdvin,
378  $ rconde, rcnde1, rcdein, scale, scale1, result,
379  $ work, lwork, iwork, info )
380 *
381 * -- LAPACK test routine (version 3.4.0) --
382 * -- LAPACK is a software package provided by Univ. of Tennessee, --
383 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
384 * November 2011
385 *
386 * .. Scalar Arguments ..
387  LOGICAL comp
388  CHARACTER balanc
389  INTEGER info, jtype, lda, ldlre, ldvl, ldvr, lwork, n,
390  $ nounit
391  DOUBLE PRECISION thresh
392 * ..
393 * .. Array Arguments ..
394  INTEGER iseed( 4 ), iwork( * )
395  DOUBLE PRECISION a( lda, * ), h( lda, * ), lre( ldlre, * ),
396  $ rcdein( * ), rcdvin( * ), rcnde1( * ),
397  $ rcndv1( * ), rconde( * ), rcondv( * ),
398  $ result( 11 ), scale( * ), scale1( * ),
399  $ vl( ldvl, * ), vr( ldvr, * ), wi( * ),
400  $ wi1( * ), work( * ), wr( * ), wr1( * )
401 * ..
402 *
403 * =====================================================================
404 *
405 *
406 * .. Parameters ..
407  DOUBLE PRECISION zero, one, two
408  parameter( zero = 0.0d0, one = 1.0d0, two = 2.0d0 )
409  DOUBLE PRECISION epsin
410  parameter( epsin = 5.9605d-8 )
411 * ..
412 * .. Local Scalars ..
413  LOGICAL balok, nobal
414  CHARACTER sense
415  INTEGER i, ihi, ihi1, iinfo, ilo, ilo1, isens, isensm,
416  $ j, jj, kmin
417  DOUBLE PRECISION abnrm, abnrm1, eps, smlnum, tnrm, tol, tolin,
418  $ ulp, ulpinv, v, vimin, vmax, vmx, vrmin, vrmx,
419  $ vtst
420 * ..
421 * .. Local Arrays ..
422  CHARACTER sens( 2 )
423  DOUBLE PRECISION dum( 1 ), res( 2 )
424 * ..
425 * .. External Functions ..
426  LOGICAL lsame
427  DOUBLE PRECISION dlamch, dlapy2, dnrm2
428  EXTERNAL lsame, dlamch, dlapy2, dnrm2
429 * ..
430 * .. External Subroutines ..
431  EXTERNAL dgeevx, dget22, dlacpy, xerbla
432 * ..
433 * .. Intrinsic Functions ..
434  INTRINSIC abs, dble, max, min
435 * ..
436 * .. Data statements ..
437  DATA sens / 'N', 'V' /
438 * ..
439 * .. Executable Statements ..
440 *
441 * Check for errors
442 *
443  nobal = lsame( balanc, 'N' )
444  balok = nobal .OR. lsame( balanc, 'P' ) .OR.
445  $ lsame( balanc, 'S' ) .OR. lsame( balanc, 'B' )
446  info = 0
447  IF( .NOT.balok ) THEN
448  info = -2
449  ELSE IF( thresh.LT.zero ) THEN
450  info = -4
451  ELSE IF( nounit.LE.0 ) THEN
452  info = -6
453  ELSE IF( n.LT.0 ) THEN
454  info = -7
455  ELSE IF( lda.LT.1 .OR. lda.LT.n ) THEN
456  info = -9
457  ELSE IF( ldvl.LT.1 .OR. ldvl.LT.n ) THEN
458  info = -16
459  ELSE IF( ldvr.LT.1 .OR. ldvr.LT.n ) THEN
460  info = -18
461  ELSE IF( ldlre.LT.1 .OR. ldlre.LT.n ) THEN
462  info = -20
463  ELSE IF( lwork.LT.3*n .OR. ( comp .AND. lwork.LT.6*n+n*n ) ) THEN
464  info = -31
465  END IF
466 *
467  IF( info.NE.0 ) THEN
468  CALL xerbla( 'DGET23', -info )
469  return
470  END IF
471 *
472 * Quick return if nothing to do
473 *
474  DO 10 i = 1, 11
475  result( i ) = -one
476  10 continue
477 *
478  IF( n.EQ.0 )
479  $ return
480 *
481 * More Important constants
482 *
483  ulp = dlamch( 'Precision' )
484  smlnum = dlamch( 'S' )
485  ulpinv = one / ulp
486 *
487 * Compute eigenvalues and eigenvectors, and test them
488 *
489  IF( lwork.GE.6*n+n*n ) THEN
490  sense = 'B'
491  isensm = 2
492  ELSE
493  sense = 'E'
494  isensm = 1
495  END IF
496  CALL dlacpy( 'F', n, n, a, lda, h, lda )
497  CALL dgeevx( balanc, 'V', 'V', sense, n, h, lda, wr, wi, vl, ldvl,
498  $ vr, ldvr, ilo, ihi, scale, abnrm, rconde, rcondv,
499  $ work, lwork, iwork, iinfo )
500  IF( iinfo.NE.0 ) THEN
501  result( 1 ) = ulpinv
502  IF( jtype.NE.22 ) THEN
503  WRITE( nounit, fmt = 9998 )'DGEEVX1', iinfo, n, jtype,
504  $ balanc, iseed
505  ELSE
506  WRITE( nounit, fmt = 9999 )'DGEEVX1', iinfo, n, iseed( 1 )
507  END IF
508  info = abs( iinfo )
509  return
510  END IF
511 *
512 * Do Test (1)
513 *
514  CALL dget22( 'N', 'N', 'N', n, a, lda, vr, ldvr, wr, wi, work,
515  $ res )
516  result( 1 ) = res( 1 )
517 *
518 * Do Test (2)
519 *
520  CALL dget22( 'T', 'N', 'T', n, a, lda, vl, ldvl, wr, wi, work,
521  $ res )
522  result( 2 ) = res( 1 )
523 *
524 * Do Test (3)
525 *
526  DO 30 j = 1, n
527  tnrm = one
528  IF( wi( j ).EQ.zero ) THEN
529  tnrm = dnrm2( n, vr( 1, j ), 1 )
530  ELSE IF( wi( j ).GT.zero ) THEN
531  tnrm = dlapy2( dnrm2( n, vr( 1, j ), 1 ),
532  $ dnrm2( n, vr( 1, j+1 ), 1 ) )
533  END IF
534  result( 3 ) = max( result( 3 ),
535  $ min( ulpinv, abs( tnrm-one ) / ulp ) )
536  IF( wi( j ).GT.zero ) THEN
537  vmx = zero
538  vrmx = zero
539  DO 20 jj = 1, n
540  vtst = dlapy2( vr( jj, j ), vr( jj, j+1 ) )
541  IF( vtst.GT.vmx )
542  $ vmx = vtst
543  IF( vr( jj, j+1 ).EQ.zero .AND. abs( vr( jj, j ) ).GT.
544  $ vrmx )vrmx = abs( vr( jj, j ) )
545  20 continue
546  IF( vrmx / vmx.LT.one-two*ulp )
547  $ result( 3 ) = ulpinv
548  END IF
549  30 continue
550 *
551 * Do Test (4)
552 *
553  DO 50 j = 1, n
554  tnrm = one
555  IF( wi( j ).EQ.zero ) THEN
556  tnrm = dnrm2( n, vl( 1, j ), 1 )
557  ELSE IF( wi( j ).GT.zero ) THEN
558  tnrm = dlapy2( dnrm2( n, vl( 1, j ), 1 ),
559  $ dnrm2( n, vl( 1, j+1 ), 1 ) )
560  END IF
561  result( 4 ) = max( result( 4 ),
562  $ min( ulpinv, abs( tnrm-one ) / ulp ) )
563  IF( wi( j ).GT.zero ) THEN
564  vmx = zero
565  vrmx = zero
566  DO 40 jj = 1, n
567  vtst = dlapy2( vl( jj, j ), vl( jj, j+1 ) )
568  IF( vtst.GT.vmx )
569  $ vmx = vtst
570  IF( vl( jj, j+1 ).EQ.zero .AND. abs( vl( jj, j ) ).GT.
571  $ vrmx )vrmx = abs( vl( jj, j ) )
572  40 continue
573  IF( vrmx / vmx.LT.one-two*ulp )
574  $ result( 4 ) = ulpinv
575  END IF
576  50 continue
577 *
578 * Test for all options of computing condition numbers
579 *
580  DO 200 isens = 1, isensm
581 *
582  sense = sens( isens )
583 *
584 * Compute eigenvalues only, and test them
585 *
586  CALL dlacpy( 'F', n, n, a, lda, h, lda )
587  CALL dgeevx( balanc, 'N', 'N', sense, n, h, lda, wr1, wi1, dum,
588  $ 1, dum, 1, ilo1, ihi1, scale1, abnrm1, rcnde1,
589  $ rcndv1, work, lwork, iwork, iinfo )
590  IF( iinfo.NE.0 ) THEN
591  result( 1 ) = ulpinv
592  IF( jtype.NE.22 ) THEN
593  WRITE( nounit, fmt = 9998 )'DGEEVX2', iinfo, n, jtype,
594  $ balanc, iseed
595  ELSE
596  WRITE( nounit, fmt = 9999 )'DGEEVX2', iinfo, n,
597  $ iseed( 1 )
598  END IF
599  info = abs( iinfo )
600  go to 190
601  END IF
602 *
603 * Do Test (5)
604 *
605  DO 60 j = 1, n
606  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
607  $ result( 5 ) = ulpinv
608  60 continue
609 *
610 * Do Test (8)
611 *
612  IF( .NOT.nobal ) THEN
613  DO 70 j = 1, n
614  IF( scale( j ).NE.scale1( j ) )
615  $ result( 8 ) = ulpinv
616  70 continue
617  IF( ilo.NE.ilo1 )
618  $ result( 8 ) = ulpinv
619  IF( ihi.NE.ihi1 )
620  $ result( 8 ) = ulpinv
621  IF( abnrm.NE.abnrm1 )
622  $ result( 8 ) = ulpinv
623  END IF
624 *
625 * Do Test (9)
626 *
627  IF( isens.EQ.2 .AND. n.GT.1 ) THEN
628  DO 80 j = 1, n
629  IF( rcondv( j ).NE.rcndv1( j ) )
630  $ result( 9 ) = ulpinv
631  80 continue
632  END IF
633 *
634 * Compute eigenvalues and right eigenvectors, and test them
635 *
636  CALL dlacpy( 'F', n, n, a, lda, h, lda )
637  CALL dgeevx( balanc, 'N', 'V', sense, n, h, lda, wr1, wi1, dum,
638  $ 1, lre, ldlre, ilo1, ihi1, scale1, abnrm1, rcnde1,
639  $ rcndv1, work, lwork, iwork, iinfo )
640  IF( iinfo.NE.0 ) THEN
641  result( 1 ) = ulpinv
642  IF( jtype.NE.22 ) THEN
643  WRITE( nounit, fmt = 9998 )'DGEEVX3', iinfo, n, jtype,
644  $ balanc, iseed
645  ELSE
646  WRITE( nounit, fmt = 9999 )'DGEEVX3', iinfo, n,
647  $ iseed( 1 )
648  END IF
649  info = abs( iinfo )
650  go to 190
651  END IF
652 *
653 * Do Test (5) again
654 *
655  DO 90 j = 1, n
656  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
657  $ result( 5 ) = ulpinv
658  90 continue
659 *
660 * Do Test (6)
661 *
662  DO 110 j = 1, n
663  DO 100 jj = 1, n
664  IF( vr( j, jj ).NE.lre( j, jj ) )
665  $ result( 6 ) = ulpinv
666  100 continue
667  110 continue
668 *
669 * Do Test (8) again
670 *
671  IF( .NOT.nobal ) THEN
672  DO 120 j = 1, n
673  IF( scale( j ).NE.scale1( j ) )
674  $ result( 8 ) = ulpinv
675  120 continue
676  IF( ilo.NE.ilo1 )
677  $ result( 8 ) = ulpinv
678  IF( ihi.NE.ihi1 )
679  $ result( 8 ) = ulpinv
680  IF( abnrm.NE.abnrm1 )
681  $ result( 8 ) = ulpinv
682  END IF
683 *
684 * Do Test (9) again
685 *
686  IF( isens.EQ.2 .AND. n.GT.1 ) THEN
687  DO 130 j = 1, n
688  IF( rcondv( j ).NE.rcndv1( j ) )
689  $ result( 9 ) = ulpinv
690  130 continue
691  END IF
692 *
693 * Compute eigenvalues and left eigenvectors, and test them
694 *
695  CALL dlacpy( 'F', n, n, a, lda, h, lda )
696  CALL dgeevx( balanc, 'V', 'N', sense, n, h, lda, wr1, wi1, lre,
697  $ ldlre, dum, 1, ilo1, ihi1, scale1, abnrm1, rcnde1,
698  $ rcndv1, work, lwork, iwork, iinfo )
699  IF( iinfo.NE.0 ) THEN
700  result( 1 ) = ulpinv
701  IF( jtype.NE.22 ) THEN
702  WRITE( nounit, fmt = 9998 )'DGEEVX4', iinfo, n, jtype,
703  $ balanc, iseed
704  ELSE
705  WRITE( nounit, fmt = 9999 )'DGEEVX4', iinfo, n,
706  $ iseed( 1 )
707  END IF
708  info = abs( iinfo )
709  go to 190
710  END IF
711 *
712 * Do Test (5) again
713 *
714  DO 140 j = 1, n
715  IF( wr( j ).NE.wr1( j ) .OR. wi( j ).NE.wi1( j ) )
716  $ result( 5 ) = ulpinv
717  140 continue
718 *
719 * Do Test (7)
720 *
721  DO 160 j = 1, n
722  DO 150 jj = 1, n
723  IF( vl( j, jj ).NE.lre( j, jj ) )
724  $ result( 7 ) = ulpinv
725  150 continue
726  160 continue
727 *
728 * Do Test (8) again
729 *
730  IF( .NOT.nobal ) THEN
731  DO 170 j = 1, n
732  IF( scale( j ).NE.scale1( j ) )
733  $ result( 8 ) = ulpinv
734  170 continue
735  IF( ilo.NE.ilo1 )
736  $ result( 8 ) = ulpinv
737  IF( ihi.NE.ihi1 )
738  $ result( 8 ) = ulpinv
739  IF( abnrm.NE.abnrm1 )
740  $ result( 8 ) = ulpinv
741  END IF
742 *
743 * Do Test (9) again
744 *
745  IF( isens.EQ.2 .AND. n.GT.1 ) THEN
746  DO 180 j = 1, n
747  IF( rcondv( j ).NE.rcndv1( j ) )
748  $ result( 9 ) = ulpinv
749  180 continue
750  END IF
751 *
752  190 continue
753 *
754  200 continue
755 *
756 * If COMP, compare condition numbers to precomputed ones
757 *
758  IF( comp ) THEN
759  CALL dlacpy( 'F', n, n, a, lda, h, lda )
760  CALL dgeevx( 'N', 'V', 'V', 'B', n, h, lda, wr, wi, vl, ldvl,
761  $ vr, ldvr, ilo, ihi, scale, abnrm, rconde, rcondv,
762  $ work, lwork, iwork, iinfo )
763  IF( iinfo.NE.0 ) THEN
764  result( 1 ) = ulpinv
765  WRITE( nounit, fmt = 9999 )'DGEEVX5', iinfo, n, iseed( 1 )
766  info = abs( iinfo )
767  go to 250
768  END IF
769 *
770 * Sort eigenvalues and condition numbers lexicographically
771 * to compare with inputs
772 *
773  DO 220 i = 1, n - 1
774  kmin = i
775  vrmin = wr( i )
776  vimin = wi( i )
777  DO 210 j = i + 1, n
778  IF( wr( j ).LT.vrmin ) THEN
779  kmin = j
780  vrmin = wr( j )
781  vimin = wi( j )
782  END IF
783  210 continue
784  wr( kmin ) = wr( i )
785  wi( kmin ) = wi( i )
786  wr( i ) = vrmin
787  wi( i ) = vimin
788  vrmin = rconde( kmin )
789  rconde( kmin ) = rconde( i )
790  rconde( i ) = vrmin
791  vrmin = rcondv( kmin )
792  rcondv( kmin ) = rcondv( i )
793  rcondv( i ) = vrmin
794  220 continue
795 *
796 * Compare condition numbers for eigenvectors
797 * taking their condition numbers into account
798 *
799  result( 10 ) = zero
800  eps = max( epsin, ulp )
801  v = max( dble( n )*eps*abnrm, smlnum )
802  IF( abnrm.EQ.zero )
803  $ v = one
804  DO 230 i = 1, n
805  IF( v.GT.rcondv( i )*rconde( i ) ) THEN
806  tol = rcondv( i )
807  ELSE
808  tol = v / rconde( i )
809  END IF
810  IF( v.GT.rcdvin( i )*rcdein( i ) ) THEN
811  tolin = rcdvin( i )
812  ELSE
813  tolin = v / rcdein( i )
814  END IF
815  tol = max( tol, smlnum / eps )
816  tolin = max( tolin, smlnum / eps )
817  IF( eps*( rcdvin( i )-tolin ).GT.rcondv( i )+tol ) THEN
818  vmax = one / eps
819  ELSE IF( rcdvin( i )-tolin.GT.rcondv( i )+tol ) THEN
820  vmax = ( rcdvin( i )-tolin ) / ( rcondv( i )+tol )
821  ELSE IF( rcdvin( i )+tolin.LT.eps*( rcondv( i )-tol ) ) THEN
822  vmax = one / eps
823  ELSE IF( rcdvin( i )+tolin.LT.rcondv( i )-tol ) THEN
824  vmax = ( rcondv( i )-tol ) / ( rcdvin( i )+tolin )
825  ELSE
826  vmax = one
827  END IF
828  result( 10 ) = max( result( 10 ), vmax )
829  230 continue
830 *
831 * Compare condition numbers for eigenvalues
832 * taking their condition numbers into account
833 *
834  result( 11 ) = zero
835  DO 240 i = 1, n
836  IF( v.GT.rcondv( i ) ) THEN
837  tol = one
838  ELSE
839  tol = v / rcondv( i )
840  END IF
841  IF( v.GT.rcdvin( i ) ) THEN
842  tolin = one
843  ELSE
844  tolin = v / rcdvin( i )
845  END IF
846  tol = max( tol, smlnum / eps )
847  tolin = max( tolin, smlnum / eps )
848  IF( eps*( rcdein( i )-tolin ).GT.rconde( i )+tol ) THEN
849  vmax = one / eps
850  ELSE IF( rcdein( i )-tolin.GT.rconde( i )+tol ) THEN
851  vmax = ( rcdein( i )-tolin ) / ( rconde( i )+tol )
852  ELSE IF( rcdein( i )+tolin.LT.eps*( rconde( i )-tol ) ) THEN
853  vmax = one / eps
854  ELSE IF( rcdein( i )+tolin.LT.rconde( i )-tol ) THEN
855  vmax = ( rconde( i )-tol ) / ( rcdein( i )+tolin )
856  ELSE
857  vmax = one
858  END IF
859  result( 11 ) = max( result( 11 ), vmax )
860  240 continue
861  250 continue
862 *
863  END IF
864 *
865  9999 format( ' DGET23: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
866  $ i6, ', INPUT EXAMPLE NUMBER = ', i4 )
867  9998 format( ' DGET23: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
868  $ i6, ', JTYPE=', i6, ', BALANC = ', a, ', ISEED=(',
869  $ 3( i5, ',' ), i5, ')' )
870 *
871  return
872 *
873 * End of DGET23
874 *
875  END