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*> \ingroup double_eig
371*
372*  =====================================================================
373      SUBROUTINE DGET23( COMP, BALANC, JTYPE, THRESH, ISEED, NOUNIT, N,
374     $                   A, LDA, H, WR, WI, WR1, WI1, VL, LDVL, VR,
375     $                   LDVR, LRE, LDLRE, RCONDV, RCNDV1, RCDVIN,
376     $                   RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1, RESULT,
377     $                   WORK, LWORK, IWORK, INFO )
378*
379*  -- LAPACK test routine --
380*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
381*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
382*
383*     .. Scalar Arguments ..
384      LOGICAL            COMP
385      CHARACTER          BALANC
386      INTEGER            INFO, JTYPE, LDA, LDLRE, LDVL, LDVR, LWORK, N,
387     $                   NOUNIT
388      DOUBLE PRECISION   THRESH
389*     ..
390*     .. Array Arguments ..
391      INTEGER            ISEED( 4 ), IWORK( * )
392      DOUBLE PRECISION   A( LDA, * ), H( LDA, * ), LRE( LDLRE, * ),
393     $                   RCDEIN( * ), RCDVIN( * ), RCNDE1( * ),
394     $                   RCNDV1( * ), RCONDE( * ), RCONDV( * ),
395     $                   RESULT( 11 ), SCALE( * ), SCALE1( * ),
396     $                   VL( LDVL, * ), VR( LDVR, * ), WI( * ),
397     $                   WI1( * ), WORK( * ), WR( * ), WR1( * )
398*     ..
399*
400*  =====================================================================
401*
402*
403*     .. Parameters ..
404      DOUBLE PRECISION   ZERO, ONE, TWO
405      PARAMETER          ( ZERO = 0.0D0, ONE = 1.0D0, TWO = 2.0D0 )
406      DOUBLE PRECISION   EPSIN
407      PARAMETER          ( EPSIN = 5.9605D-8 )
408*     ..
409*     .. Local Scalars ..
410      LOGICAL            BALOK, NOBAL
411      CHARACTER          SENSE
412      INTEGER            I, IHI, IHI1, IINFO, ILO, ILO1, ISENS, ISENSM,
413     $                   J, JJ, KMIN
414      DOUBLE PRECISION   ABNRM, ABNRM1, EPS, SMLNUM, TNRM, TOL, TOLIN,
415     $                   ULP, ULPINV, V, VIMIN, VMAX, VMX, VRMIN, VRMX,
416     $                   VTST
417*     ..
418*     .. Local Arrays ..
419      CHARACTER          SENS( 2 )
420      DOUBLE PRECISION   DUM( 1 ), RES( 2 )
421*     ..
422*     .. External Functions ..
423      LOGICAL            LSAME
424      DOUBLE PRECISION   DLAMCH, DLAPY2, DNRM2
425      EXTERNAL           LSAME, DLAMCH, DLAPY2, DNRM2
426*     ..
427*     .. External Subroutines ..
428      EXTERNAL           DGEEVX, DGET22, DLACPY, XERBLA
429*     ..
430*     .. Intrinsic Functions ..
431      INTRINSIC          ABS, DBLE, MAX, MIN
432*     ..
433*     .. Data statements ..
434      DATA               SENS / 'N', 'V' /
435*     ..
436*     .. Executable Statements ..
437*
438*     Check for errors
439*
440      NOBAL = LSAME( BALANC, 'N' )
441      BALOK = NOBAL .OR. LSAME( BALANC, 'P' ) .OR.
442     $        LSAME( BALANC, 'S' ) .OR. LSAME( BALANC, 'B' )
443      INFO = 0
444      IF( .NOT.BALOK ) THEN
445         INFO = -2
446      ELSE IF( THRESH.LT.ZERO ) THEN
447         INFO = -4
448      ELSE IF( NOUNIT.LE.0 ) THEN
449         INFO = -6
450      ELSE IF( N.LT.0 ) THEN
451         INFO = -7
452      ELSE IF( LDA.LT.1 .OR. LDA.LT.N ) THEN
453         INFO = -9
454      ELSE IF( LDVL.LT.1 .OR. LDVL.LT.N ) THEN
455         INFO = -16
456      ELSE IF( LDVR.LT.1 .OR. LDVR.LT.N ) THEN
457         INFO = -18
458      ELSE IF( LDLRE.LT.1 .OR. LDLRE.LT.N ) THEN
459         INFO = -20
460      ELSE IF( LWORK.LT.3*N .OR. ( COMP .AND. LWORK.LT.6*N+N*N ) ) THEN
461         INFO = -31
462      END IF
463*
464      IF( INFO.NE.0 ) THEN
465         CALL XERBLA( 'DGET23', -INFO )
466         RETURN
467      END IF
468*
469*     Quick return if nothing to do
470*
471      DO 10 I = 1, 11
472         RESULT( I ) = -ONE
473   10 CONTINUE
474*
475      IF( N.EQ.0 )
476     $   RETURN
477*
478*     More Important constants
479*
480      ULP = DLAMCH( 'Precision' )
481      SMLNUM = DLAMCH( 'S' )
482      ULPINV = ONE / ULP
483*
484*     Compute eigenvalues and eigenvectors, and test them
485*
486      IF( LWORK.GE.6*N+N*N ) THEN
487         SENSE = 'B'
488         ISENSM = 2
489      ELSE
490         SENSE = 'E'
491         ISENSM = 1
492      END IF
493      CALL DLACPY( 'F', N, N, A, LDA, H, LDA )
494      CALL DGEEVX( BALANC, 'V', 'V', SENSE, N, H, LDA, WR, WI, VL, LDVL,
495     $             VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE, RCONDV,
496     $             WORK, LWORK, IWORK, IINFO )
497      IF( IINFO.NE.0 ) THEN
498         RESULT( 1 ) = ULPINV
499         IF( JTYPE.NE.22 ) THEN
500            WRITE( NOUNIT, FMT = 9998 )'DGEEVX1', IINFO, N, JTYPE,
501     $         BALANC, ISEED
502         ELSE
503            WRITE( NOUNIT, FMT = 9999 )'DGEEVX1', IINFO, N, ISEED( 1 )
504         END IF
505         INFO = ABS( IINFO )
506         RETURN
507      END IF
508*
509*     Do Test (1)
510*
511      CALL DGET22( 'N', 'N', 'N', N, A, LDA, VR, LDVR, WR, WI, WORK,
512     $             RES )
513      RESULT( 1 ) = RES( 1 )
514*
515*     Do Test (2)
516*
517      CALL DGET22( 'T', 'N', 'T', N, A, LDA, VL, LDVL, WR, WI, WORK,
518     $             RES )
519      RESULT( 2 ) = RES( 1 )
520*
521*     Do Test (3)
522*
523      DO 30 J = 1, N
524         TNRM = ONE
525         IF( WI( J ).EQ.ZERO ) THEN
526            TNRM = DNRM2( N, VR( 1, J ), 1 )
527         ELSE IF( WI( J ).GT.ZERO ) THEN
528            TNRM = DLAPY2( DNRM2( N, VR( 1, J ), 1 ),
529     $             DNRM2( N, VR( 1, J+1 ), 1 ) )
530         END IF
531         RESULT( 3 ) = MAX( RESULT( 3 ),
532     $                 MIN( ULPINV, ABS( TNRM-ONE ) / ULP ) )
533         IF( WI( J ).GT.ZERO ) THEN
534            VMX = ZERO
535            VRMX = ZERO
536            DO 20 JJ = 1, N
537               VTST = DLAPY2( VR( JJ, J ), VR( JJ, J+1 ) )
538               IF( VTST.GT.VMX )
539     $            VMX = VTST
540               IF( VR( JJ, J+1 ).EQ.ZERO .AND. ABS( VR( JJ, J ) ).GT.
541     $             VRMX )VRMX = ABS( VR( JJ, J ) )
542   20       CONTINUE
543            IF( VRMX / VMX.LT.ONE-TWO*ULP )
544     $         RESULT( 3 ) = ULPINV
545         END IF
546   30 CONTINUE
547*
548*     Do Test (4)
549*
550      DO 50 J = 1, N
551         TNRM = ONE
552         IF( WI( J ).EQ.ZERO ) THEN
553            TNRM = DNRM2( N, VL( 1, J ), 1 )
554         ELSE IF( WI( J ).GT.ZERO ) THEN
555            TNRM = DLAPY2( DNRM2( N, VL( 1, J ), 1 ),
556     $             DNRM2( N, VL( 1, J+1 ), 1 ) )
557         END IF
558         RESULT( 4 ) = MAX( RESULT( 4 ),
559     $                 MIN( ULPINV, ABS( TNRM-ONE ) / ULP ) )
560         IF( WI( J ).GT.ZERO ) THEN
561            VMX = ZERO
562            VRMX = ZERO
563            DO 40 JJ = 1, N
564               VTST = DLAPY2( VL( JJ, J ), VL( JJ, J+1 ) )
565               IF( VTST.GT.VMX )
566     $            VMX = VTST
567               IF( VL( JJ, J+1 ).EQ.ZERO .AND. ABS( VL( JJ, J ) ).GT.
568     $             VRMX )VRMX = ABS( VL( JJ, J ) )
569   40       CONTINUE
570            IF( VRMX / VMX.LT.ONE-TWO*ULP )
571     $         RESULT( 4 ) = ULPINV
572         END IF
573   50 CONTINUE
574*
575*     Test for all options of computing condition numbers
576*
577      DO 200 ISENS = 1, ISENSM
578*
579         SENSE = SENS( ISENS )
580*
581*        Compute eigenvalues only, and test them
582*
583         CALL DLACPY( 'F', N, N, A, LDA, H, LDA )
584         CALL DGEEVX( BALANC, 'N', 'N', SENSE, N, H, LDA, WR1, WI1, DUM,
585     $                1, DUM, 1, ILO1, IHI1, SCALE1, ABNRM1, RCNDE1,
586     $                RCNDV1, WORK, LWORK, IWORK, IINFO )
587         IF( IINFO.NE.0 ) THEN
588            RESULT( 1 ) = ULPINV
589            IF( JTYPE.NE.22 ) THEN
590               WRITE( NOUNIT, FMT = 9998 )'DGEEVX2', IINFO, N, JTYPE,
591     $            BALANC, ISEED
592            ELSE
593               WRITE( NOUNIT, FMT = 9999 )'DGEEVX2', IINFO, N,
594     $            ISEED( 1 )
595            END IF
596            INFO = ABS( IINFO )
597            GO TO 190
598         END IF
599*
600*        Do Test (5)
601*
602         DO 60 J = 1, N
603            IF( WR( J ).NE.WR1( J ) .OR. WI( J ).NE.WI1( J ) )
604     $         RESULT( 5 ) = ULPINV
605   60    CONTINUE
606*
607*        Do Test (8)
608*
609         IF( .NOT.NOBAL ) THEN
610            DO 70 J = 1, N
611               IF( SCALE( J ).NE.SCALE1( J ) )
612     $            RESULT( 8 ) = ULPINV
613   70       CONTINUE
614            IF( ILO.NE.ILO1 )
615     $         RESULT( 8 ) = ULPINV
616            IF( IHI.NE.IHI1 )
617     $         RESULT( 8 ) = ULPINV
618            IF( ABNRM.NE.ABNRM1 )
619     $         RESULT( 8 ) = ULPINV
620         END IF
621*
622*        Do Test (9)
623*
624         IF( ISENS.EQ.2 .AND. N.GT.1 ) THEN
625            DO 80 J = 1, N
626               IF( RCONDV( J ).NE.RCNDV1( J ) )
627     $            RESULT( 9 ) = ULPINV
628   80       CONTINUE
629         END IF
630*
631*        Compute eigenvalues and right eigenvectors, and test them
632*
633         CALL DLACPY( 'F', N, N, A, LDA, H, LDA )
634         CALL DGEEVX( BALANC, 'N', 'V', SENSE, N, H, LDA, WR1, WI1, DUM,
635     $                1, LRE, LDLRE, ILO1, IHI1, SCALE1, ABNRM1, RCNDE1,
636     $                RCNDV1, WORK, LWORK, IWORK, IINFO )
637         IF( IINFO.NE.0 ) THEN
638            RESULT( 1 ) = ULPINV
639            IF( JTYPE.NE.22 ) THEN
640               WRITE( NOUNIT, FMT = 9998 )'DGEEVX3', IINFO, N, JTYPE,
641     $            BALANC, ISEED
642            ELSE
643               WRITE( NOUNIT, FMT = 9999 )'DGEEVX3', IINFO, N,
644     $            ISEED( 1 )
645            END IF
646            INFO = ABS( IINFO )
647            GO TO 190
648         END IF
649*
650*        Do Test (5) again
651*
652         DO 90 J = 1, N
653            IF( WR( J ).NE.WR1( J ) .OR. WI( J ).NE.WI1( J ) )
654     $         RESULT( 5 ) = ULPINV
655   90    CONTINUE
656*
657*        Do Test (6)
658*
659         DO 110 J = 1, N
660            DO 100 JJ = 1, N
661               IF( VR( J, JJ ).NE.LRE( J, JJ ) )
662     $            RESULT( 6 ) = ULPINV
663  100       CONTINUE
664  110    CONTINUE
665*
666*        Do Test (8) again
667*
668         IF( .NOT.NOBAL ) THEN
669            DO 120 J = 1, N
670               IF( SCALE( J ).NE.SCALE1( J ) )
671     $            RESULT( 8 ) = ULPINV
672  120       CONTINUE
673            IF( ILO.NE.ILO1 )
674     $         RESULT( 8 ) = ULPINV
675            IF( IHI.NE.IHI1 )
676     $         RESULT( 8 ) = ULPINV
677            IF( ABNRM.NE.ABNRM1 )
678     $         RESULT( 8 ) = ULPINV
679         END IF
680*
681*        Do Test (9) again
682*
683         IF( ISENS.EQ.2 .AND. N.GT.1 ) THEN
684            DO 130 J = 1, N
685               IF( RCONDV( J ).NE.RCNDV1( J ) )
686     $            RESULT( 9 ) = ULPINV
687  130       CONTINUE
688         END IF
689*
690*        Compute eigenvalues and left eigenvectors, and test them
691*
692         CALL DLACPY( 'F', N, N, A, LDA, H, LDA )
693         CALL DGEEVX( BALANC, 'V', 'N', SENSE, N, H, LDA, WR1, WI1, LRE,
694     $                LDLRE, DUM, 1, ILO1, IHI1, SCALE1, ABNRM1, RCNDE1,
695     $                RCNDV1, WORK, LWORK, IWORK, IINFO )
696         IF( IINFO.NE.0 ) THEN
697            RESULT( 1 ) = ULPINV
698            IF( JTYPE.NE.22 ) THEN
699               WRITE( NOUNIT, FMT = 9998 )'DGEEVX4', IINFO, N, JTYPE,
700     $            BALANC, ISEED
701            ELSE
702               WRITE( NOUNIT, FMT = 9999 )'DGEEVX4', IINFO, N,
703     $            ISEED( 1 )
704            END IF
705            INFO = ABS( IINFO )
706            GO TO 190
707         END IF
708*
709*        Do Test (5) again
710*
711         DO 140 J = 1, N
712            IF( WR( J ).NE.WR1( J ) .OR. WI( J ).NE.WI1( J ) )
713     $         RESULT( 5 ) = ULPINV
714  140    CONTINUE
715*
716*        Do Test (7)
717*
718         DO 160 J = 1, N
719            DO 150 JJ = 1, N
720               IF( VL( J, JJ ).NE.LRE( J, JJ ) )
721     $            RESULT( 7 ) = ULPINV
722  150       CONTINUE
723  160    CONTINUE
724*
725*        Do Test (8) again
726*
727         IF( .NOT.NOBAL ) THEN
728            DO 170 J = 1, N
729               IF( SCALE( J ).NE.SCALE1( J ) )
730     $            RESULT( 8 ) = ULPINV
731  170       CONTINUE
732            IF( ILO.NE.ILO1 )
733     $         RESULT( 8 ) = ULPINV
734            IF( IHI.NE.IHI1 )
735     $         RESULT( 8 ) = ULPINV
736            IF( ABNRM.NE.ABNRM1 )
737     $         RESULT( 8 ) = ULPINV
738         END IF
739*
740*        Do Test (9) again
741*
742         IF( ISENS.EQ.2 .AND. N.GT.1 ) THEN
743            DO 180 J = 1, N
744               IF( RCONDV( J ).NE.RCNDV1( J ) )
745     $            RESULT( 9 ) = ULPINV
746  180       CONTINUE
747         END IF
748*
749  190    CONTINUE
750*
751  200 CONTINUE
752*
753*     If COMP, compare condition numbers to precomputed ones
754*
755      IF( COMP ) THEN
756         CALL DLACPY( 'F', N, N, A, LDA, H, LDA )
757         CALL DGEEVX( 'N', 'V', 'V', 'B', N, H, LDA, WR, WI, VL, LDVL,
758     $                VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE, RCONDV,
759     $                WORK, LWORK, IWORK, IINFO )
760         IF( IINFO.NE.0 ) THEN
761            RESULT( 1 ) = ULPINV
762            WRITE( NOUNIT, FMT = 9999 )'DGEEVX5', IINFO, N, ISEED( 1 )
763            INFO = ABS( IINFO )
764            GO TO 250
765         END IF
766*
767*        Sort eigenvalues and condition numbers lexicographically
768*        to compare with inputs
769*
770         DO 220 I = 1, N - 1
771            KMIN = I
772            VRMIN = WR( I )
773            VIMIN = WI( I )
774            DO 210 J = I + 1, N
775               IF( WR( J ).LT.VRMIN ) THEN
776                  KMIN = J
777                  VRMIN = WR( J )
778                  VIMIN = WI( J )
779               END IF
780  210       CONTINUE
781            WR( KMIN ) = WR( I )
782            WI( KMIN ) = WI( I )
783            WR( I ) = VRMIN
784            WI( I ) = VIMIN
785            VRMIN = RCONDE( KMIN )
786            RCONDE( KMIN ) = RCONDE( I )
787            RCONDE( I ) = VRMIN
788            VRMIN = RCONDV( KMIN )
789            RCONDV( KMIN ) = RCONDV( I )
790            RCONDV( I ) = VRMIN
791  220    CONTINUE
792*
793*        Compare condition numbers for eigenvectors
794*        taking their condition numbers into account
795*
796         RESULT( 10 ) = ZERO
797         EPS = MAX( EPSIN, ULP )
798         V = MAX( DBLE( N )*EPS*ABNRM, SMLNUM )
799         IF( ABNRM.EQ.ZERO )
800     $      V = ONE
801         DO 230 I = 1, N
802            IF( V.GT.RCONDV( I )*RCONDE( I ) ) THEN
803               TOL = RCONDV( I )
804            ELSE
805               TOL = V / RCONDE( I )
806            END IF
807            IF( V.GT.RCDVIN( I )*RCDEIN( I ) ) THEN
808               TOLIN = RCDVIN( I )
809            ELSE
810               TOLIN = V / RCDEIN( I )
811            END IF
812            TOL = MAX( TOL, SMLNUM / EPS )
813            TOLIN = MAX( TOLIN, SMLNUM / EPS )
814            IF( EPS*( RCDVIN( I )-TOLIN ).GT.RCONDV( I )+TOL ) THEN
815               VMAX = ONE / EPS
816            ELSE IF( RCDVIN( I )-TOLIN.GT.RCONDV( I )+TOL ) THEN
817               VMAX = ( RCDVIN( I )-TOLIN ) / ( RCONDV( I )+TOL )
818            ELSE IF( RCDVIN( I )+TOLIN.LT.EPS*( RCONDV( I )-TOL ) ) THEN
819               VMAX = ONE / EPS
820            ELSE IF( RCDVIN( I )+TOLIN.LT.RCONDV( I )-TOL ) THEN
821               VMAX = ( RCONDV( I )-TOL ) / ( RCDVIN( I )+TOLIN )
822            ELSE
823               VMAX = ONE
824            END IF
825            RESULT( 10 ) = MAX( RESULT( 10 ), VMAX )
826  230    CONTINUE
827*
828*        Compare condition numbers for eigenvalues
829*        taking their condition numbers into account
830*
831         RESULT( 11 ) = ZERO
832         DO 240 I = 1, N
833            IF( V.GT.RCONDV( I ) ) THEN
834               TOL = ONE
835            ELSE
836               TOL = V / RCONDV( I )
837            END IF
838            IF( V.GT.RCDVIN( I ) ) THEN
839               TOLIN = ONE
840            ELSE
841               TOLIN = V / RCDVIN( I )
842            END IF
843            TOL = MAX( TOL, SMLNUM / EPS )
844            TOLIN = MAX( TOLIN, SMLNUM / EPS )
845            IF( EPS*( RCDEIN( I )-TOLIN ).GT.RCONDE( I )+TOL ) THEN
846               VMAX = ONE / EPS
847            ELSE IF( RCDEIN( I )-TOLIN.GT.RCONDE( I )+TOL ) THEN
848               VMAX = ( RCDEIN( I )-TOLIN ) / ( RCONDE( I )+TOL )
849            ELSE IF( RCDEIN( I )+TOLIN.LT.EPS*( RCONDE( I )-TOL ) ) THEN
850               VMAX = ONE / EPS
851            ELSE IF( RCDEIN( I )+TOLIN.LT.RCONDE( I )-TOL ) THEN
852               VMAX = ( RCONDE( I )-TOL ) / ( RCDEIN( I )+TOLIN )
853            ELSE
854               VMAX = ONE
855            END IF
856            RESULT( 11 ) = MAX( RESULT( 11 ), VMAX )
857  240    CONTINUE
858  250    CONTINUE
859*
860      END IF
861*
862 9999 FORMAT( ' DGET23: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
863     $      I6, ', INPUT EXAMPLE NUMBER = ', I4 )
864 9998 FORMAT( ' DGET23: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
865     $      I6, ', JTYPE=', I6, ', BALANC = ', A, ', ISEED=(',
866     $      3( I5, ',' ), I5, ')' )
867*
868      RETURN
869*
870*     End of DGET23
871*
872      END
873