LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zgeevx.f
Go to the documentation of this file.
1 *> \brief <b> ZGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZGEEVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgeevx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgeevx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgeevx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZGEEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, W, VL,
22 * LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE,
23 * RCONDV, WORK, LWORK, RWORK, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER BALANC, JOBVL, JOBVR, SENSE
27 * INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
28 * DOUBLE PRECISION ABNRM
29 * ..
30 * .. Array Arguments ..
31 * DOUBLE PRECISION RCONDE( * ), RCONDV( * ), RWORK( * ),
32 * $ SCALE( * )
33 * COMPLEX*16 A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ),
34 * $ W( * ), WORK( * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> ZGEEVX computes for an N-by-N complex nonsymmetric matrix A, the
44 *> eigenvalues and, optionally, the left and/or right eigenvectors.
45 *>
46 *> Optionally also, it computes a balancing transformation to improve
47 *> the conditioning of the eigenvalues and eigenvectors (ILO, IHI,
48 *> SCALE, and ABNRM), reciprocal condition numbers for the eigenvalues
49 *> (RCONDE), and reciprocal condition numbers for the right
50 *> eigenvectors (RCONDV).
51 *>
52 *> The right eigenvector v(j) of A satisfies
53 *> A * v(j) = lambda(j) * v(j)
54 *> where lambda(j) is its eigenvalue.
55 *> The left eigenvector u(j) of A satisfies
56 *> u(j)**H * A = lambda(j) * u(j)**H
57 *> where u(j)**H denotes the conjugate transpose of u(j).
58 *>
59 *> The computed eigenvectors are normalized to have Euclidean norm
60 *> equal to 1 and largest component real.
61 *>
62 *> Balancing a matrix means permuting the rows and columns to make it
63 *> more nearly upper triangular, and applying a diagonal similarity
64 *> transformation D * A * D**(-1), where D is a diagonal matrix, to
65 *> make its rows and columns closer in norm and the condition numbers
66 *> of its eigenvalues and eigenvectors smaller. The computed
67 *> reciprocal condition numbers correspond to the balanced matrix.
68 *> Permuting rows and columns will not change the condition numbers
69 *> (in exact arithmetic) but diagonal scaling will. For further
70 *> explanation of balancing, see section 4.10.2 of the LAPACK
71 *> Users' Guide.
72 *> \endverbatim
73 *
74 * Arguments:
75 * ==========
76 *
77 *> \param[in] BALANC
78 *> \verbatim
79 *> BALANC is CHARACTER*1
80 *> Indicates how the input matrix should be diagonally scaled
81 *> and/or permuted to improve the conditioning of its
82 *> eigenvalues.
83 *> = 'N': Do not diagonally scale or permute;
84 *> = 'P': Perform permutations to make the matrix more nearly
85 *> upper triangular. Do not diagonally scale;
86 *> = 'S': Diagonally scale the matrix, ie. replace A by
87 *> D*A*D**(-1), where D is a diagonal matrix chosen
88 *> to make the rows and columns of A more equal in
89 *> norm. Do not permute;
90 *> = 'B': Both diagonally scale and permute A.
91 *>
92 *> Computed reciprocal condition numbers will be for the matrix
93 *> after balancing and/or permuting. Permuting does not change
94 *> condition numbers (in exact arithmetic), but balancing does.
95 *> \endverbatim
96 *>
97 *> \param[in] JOBVL
98 *> \verbatim
99 *> JOBVL is CHARACTER*1
100 *> = 'N': left eigenvectors of A are not computed;
101 *> = 'V': left eigenvectors of A are computed.
102 *> If SENSE = 'E' or 'B', JOBVL must = 'V'.
103 *> \endverbatim
104 *>
105 *> \param[in] JOBVR
106 *> \verbatim
107 *> JOBVR is CHARACTER*1
108 *> = 'N': right eigenvectors of A are not computed;
109 *> = 'V': right eigenvectors of A are computed.
110 *> If SENSE = 'E' or 'B', JOBVR must = 'V'.
111 *> \endverbatim
112 *>
113 *> \param[in] SENSE
114 *> \verbatim
115 *> SENSE is CHARACTER*1
116 *> Determines which reciprocal condition numbers are computed.
117 *> = 'N': None are computed;
118 *> = 'E': Computed for eigenvalues only;
119 *> = 'V': Computed for right eigenvectors only;
120 *> = 'B': Computed for eigenvalues and right eigenvectors.
121 *>
122 *> If SENSE = 'E' or 'B', both left and right eigenvectors
123 *> must also be computed (JOBVL = 'V' and JOBVR = 'V').
124 *> \endverbatim
125 *>
126 *> \param[in] N
127 *> \verbatim
128 *> N is INTEGER
129 *> The order of the matrix A. N >= 0.
130 *> \endverbatim
131 *>
132 *> \param[in,out] A
133 *> \verbatim
134 *> A is COMPLEX*16 array, dimension (LDA,N)
135 *> On entry, the N-by-N matrix A.
136 *> On exit, A has been overwritten. If JOBVL = 'V' or
137 *> JOBVR = 'V', A contains the Schur form of the balanced
138 *> version of the matrix A.
139 *> \endverbatim
140 *>
141 *> \param[in] LDA
142 *> \verbatim
143 *> LDA is INTEGER
144 *> The leading dimension of the array A. LDA >= max(1,N).
145 *> \endverbatim
146 *>
147 *> \param[out] W
148 *> \verbatim
149 *> W is COMPLEX*16 array, dimension (N)
150 *> W contains the computed eigenvalues.
151 *> \endverbatim
152 *>
153 *> \param[out] VL
154 *> \verbatim
155 *> VL is COMPLEX*16 array, dimension (LDVL,N)
156 *> If JOBVL = 'V', the left eigenvectors u(j) are stored one
157 *> after another in the columns of VL, in the same order
158 *> as their eigenvalues.
159 *> If JOBVL = 'N', VL is not referenced.
160 *> u(j) = VL(:,j), the j-th column of VL.
161 *> \endverbatim
162 *>
163 *> \param[in] LDVL
164 *> \verbatim
165 *> LDVL is INTEGER
166 *> The leading dimension of the array VL. LDVL >= 1; if
167 *> JOBVL = 'V', LDVL >= N.
168 *> \endverbatim
169 *>
170 *> \param[out] VR
171 *> \verbatim
172 *> VR is COMPLEX*16 array, dimension (LDVR,N)
173 *> If JOBVR = 'V', the right eigenvectors v(j) are stored one
174 *> after another in the columns of VR, in the same order
175 *> as their eigenvalues.
176 *> If JOBVR = 'N', VR is not referenced.
177 *> v(j) = VR(:,j), the j-th column of VR.
178 *> \endverbatim
179 *>
180 *> \param[in] LDVR
181 *> \verbatim
182 *> LDVR is INTEGER
183 *> The leading dimension of the array VR. LDVR >= 1; if
184 *> JOBVR = 'V', LDVR >= N.
185 *> \endverbatim
186 *>
187 *> \param[out] ILO
188 *> \verbatim
189 *> ILO is INTEGER
190 *> \endverbatim
191 *>
192 *> \param[out] IHI
193 *> \verbatim
194 *> IHI is INTEGER
195 *> ILO and IHI are integer values determined when A was
196 *> balanced. The balanced A(i,j) = 0 if I > J and
197 *> J = 1,...,ILO-1 or I = IHI+1,...,N.
198 *> \endverbatim
199 *>
200 *> \param[out] SCALE
201 *> \verbatim
202 *> SCALE is DOUBLE PRECISION array, dimension (N)
203 *> Details of the permutations and scaling factors applied
204 *> when balancing A. If P(j) is the index of the row and column
205 *> interchanged with row and column j, and D(j) is the scaling
206 *> factor applied to row and column j, then
207 *> SCALE(J) = P(J), for J = 1,...,ILO-1
208 *> = D(J), for J = ILO,...,IHI
209 *> = P(J) for J = IHI+1,...,N.
210 *> The order in which the interchanges are made is N to IHI+1,
211 *> then 1 to ILO-1.
212 *> \endverbatim
213 *>
214 *> \param[out] ABNRM
215 *> \verbatim
216 *> ABNRM is DOUBLE PRECISION
217 *> The one-norm of the balanced matrix (the maximum
218 *> of the sum of absolute values of elements of any column).
219 *> \endverbatim
220 *>
221 *> \param[out] RCONDE
222 *> \verbatim
223 *> RCONDE is DOUBLE PRECISION array, dimension (N)
224 *> RCONDE(j) is the reciprocal condition number of the j-th
225 *> eigenvalue.
226 *> \endverbatim
227 *>
228 *> \param[out] RCONDV
229 *> \verbatim
230 *> RCONDV is DOUBLE PRECISION array, dimension (N)
231 *> RCONDV(j) is the reciprocal condition number of the j-th
232 *> right eigenvector.
233 *> \endverbatim
234 *>
235 *> \param[out] WORK
236 *> \verbatim
237 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
238 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
239 *> \endverbatim
240 *>
241 *> \param[in] LWORK
242 *> \verbatim
243 *> LWORK is INTEGER
244 *> The dimension of the array WORK. If SENSE = 'N' or 'E',
245 *> LWORK >= max(1,2*N), and if SENSE = 'V' or 'B',
246 *> LWORK >= N*N+2*N.
247 *> For good performance, LWORK must generally be larger.
248 *>
249 *> If LWORK = -1, then a workspace query is assumed; the routine
250 *> only calculates the optimal size of the WORK array, returns
251 *> this value as the first entry of the WORK array, and no error
252 *> message related to LWORK is issued by XERBLA.
253 *> \endverbatim
254 *>
255 *> \param[out] RWORK
256 *> \verbatim
257 *> RWORK is DOUBLE PRECISION array, dimension (2*N)
258 *> \endverbatim
259 *>
260 *> \param[out] INFO
261 *> \verbatim
262 *> INFO is INTEGER
263 *> = 0: successful exit
264 *> < 0: if INFO = -i, the i-th argument had an illegal value.
265 *> > 0: if INFO = i, the QR algorithm failed to compute all the
266 *> eigenvalues, and no eigenvectors or condition numbers
267 *> have been computed; elements 1:ILO-1 and i+1:N of W
268 *> contain eigenvalues which have converged.
269 *> \endverbatim
270 *
271 * Authors:
272 * ========
273 *
274 *> \author Univ. of Tennessee
275 *> \author Univ. of California Berkeley
276 *> \author Univ. of Colorado Denver
277 *> \author NAG Ltd.
278 *
279 *> \date June 2016
280 *
281 * @precisions fortran z -> c
282 *
283 *> \ingroup complex16GEeigen
284 *
285 * =====================================================================
286  SUBROUTINE zgeevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, W, VL,
287  $ ldvl, vr, ldvr, ilo, ihi, scale, abnrm, rconde,
288  $ rcondv, work, lwork, rwork, info )
289  implicit none
290 *
291 * -- LAPACK driver routine (version 3.6.1) --
292 * -- LAPACK is a software package provided by Univ. of Tennessee, --
293 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
294 * June 2016
295 *
296 * .. Scalar Arguments ..
297  CHARACTER BALANC, JOBVL, JOBVR, SENSE
298  INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
299  DOUBLE PRECISION ABNRM
300 * ..
301 * .. Array Arguments ..
302  DOUBLE PRECISION RCONDE( * ), RCONDV( * ), RWORK( * ),
303  $ scale( * )
304  COMPLEX*16 A( lda, * ), VL( ldvl, * ), VR( ldvr, * ),
305  $ w( * ), work( * )
306 * ..
307 *
308 * =====================================================================
309 *
310 * .. Parameters ..
311  DOUBLE PRECISION ZERO, ONE
312  parameter ( zero = 0.0d0, one = 1.0d0 )
313 * ..
314 * .. Local Scalars ..
315  LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
316  $ wntsnn, wntsnv
317  CHARACTER JOB, SIDE
318  INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
319  $ lwork_trevc, maxwrk, minwrk, nout
320  DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, SCL, SMLNUM
321  COMPLEX*16 TMP
322 * ..
323 * .. Local Arrays ..
324  LOGICAL SELECT( 1 )
325  DOUBLE PRECISION DUM( 1 )
326 * ..
327 * .. External Subroutines ..
328  EXTERNAL dlabad, dlascl, xerbla, zdscal, zgebak, zgebal,
330  $ ztrsna, zunghr
331 * ..
332 * .. External Functions ..
333  LOGICAL LSAME
334  INTEGER IDAMAX, ILAENV
335  DOUBLE PRECISION DLAMCH, DZNRM2, ZLANGE
336  EXTERNAL lsame, idamax, ilaenv, dlamch, dznrm2, zlange
337 * ..
338 * .. Intrinsic Functions ..
339  INTRINSIC dble, dcmplx, conjg, aimag, max, sqrt
340 * ..
341 * .. Executable Statements ..
342 *
343 * Test the input arguments
344 *
345  info = 0
346  lquery = ( lwork.EQ.-1 )
347  wantvl = lsame( jobvl, 'V' )
348  wantvr = lsame( jobvr, 'V' )
349  wntsnn = lsame( sense, 'N' )
350  wntsne = lsame( sense, 'E' )
351  wntsnv = lsame( sense, 'V' )
352  wntsnb = lsame( sense, 'B' )
353  IF( .NOT.( lsame( balanc, 'N' ) .OR. lsame( balanc, 'S' ) .OR.
354  $ lsame( balanc, 'P' ) .OR. lsame( balanc, 'B' ) ) ) THEN
355  info = -1
356  ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
357  info = -2
358  ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
359  info = -3
360  ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
361  $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
362  $ wantvr ) ) ) THEN
363  info = -4
364  ELSE IF( n.LT.0 ) THEN
365  info = -5
366  ELSE IF( lda.LT.max( 1, n ) ) THEN
367  info = -7
368  ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
369  info = -10
370  ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
371  info = -12
372  END IF
373 *
374 * Compute workspace
375 * (Note: Comments in the code beginning "Workspace:" describe the
376 * minimal amount of workspace needed at that point in the code,
377 * as well as the preferred amount for good performance.
378 * CWorkspace refers to complex workspace, and RWorkspace to real
379 * workspace. NB refers to the optimal block size for the
380 * immediately following subroutine, as returned by ILAENV.
381 * HSWORK refers to the workspace preferred by ZHSEQR, as
382 * calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
383 * the worst case.)
384 *
385  IF( info.EQ.0 ) THEN
386  IF( n.EQ.0 ) THEN
387  minwrk = 1
388  maxwrk = 1
389  ELSE
390  maxwrk = n + n*ilaenv( 1, 'ZGEHRD', ' ', n, 1, n, 0 )
391 *
392  IF( wantvl ) THEN
393  CALL ztrevc3( 'L', 'B', SELECT, n, a, lda,
394  $ vl, ldvl, vr, ldvr,
395  $ n, nout, work, -1, rwork, -1, ierr )
396  lwork_trevc = int( work(1) )
397  maxwrk = max( maxwrk, lwork_trevc )
398  CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vl, ldvl,
399  $ work, -1, info )
400  ELSE IF( wantvr ) THEN
401  CALL ztrevc3( 'R', 'B', SELECT, n, a, lda,
402  $ vl, ldvl, vr, ldvr,
403  $ n, nout, work, -1, rwork, -1, ierr )
404  lwork_trevc = int( work(1) )
405  maxwrk = max( maxwrk, lwork_trevc )
406  CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vr, ldvr,
407  $ work, -1, info )
408  ELSE
409  IF( wntsnn ) THEN
410  CALL zhseqr( 'E', 'N', n, 1, n, a, lda, w, vr, ldvr,
411  $ work, -1, info )
412  ELSE
413  CALL zhseqr( 'S', 'N', n, 1, n, a, lda, w, vr, ldvr,
414  $ work, -1, info )
415  END IF
416  END IF
417  hswork = int( work(1) )
418 *
419  IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) ) THEN
420  minwrk = 2*n
421  IF( .NOT.( wntsnn .OR. wntsne ) )
422  $ minwrk = max( minwrk, n*n + 2*n )
423  maxwrk = max( maxwrk, hswork )
424  IF( .NOT.( wntsnn .OR. wntsne ) )
425  $ maxwrk = max( maxwrk, n*n + 2*n )
426  ELSE
427  minwrk = 2*n
428  IF( .NOT.( wntsnn .OR. wntsne ) )
429  $ minwrk = max( minwrk, n*n + 2*n )
430  maxwrk = max( maxwrk, hswork )
431  maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1, 'ZUNGHR',
432  $ ' ', n, 1, n, -1 ) )
433  IF( .NOT.( wntsnn .OR. wntsne ) )
434  $ maxwrk = max( maxwrk, n*n + 2*n )
435  maxwrk = max( maxwrk, 2*n )
436  END IF
437  maxwrk = max( maxwrk, minwrk )
438  END IF
439  work( 1 ) = maxwrk
440 *
441  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
442  info = -20
443  END IF
444  END IF
445 *
446  IF( info.NE.0 ) THEN
447  CALL xerbla( 'ZGEEVX', -info )
448  RETURN
449  ELSE IF( lquery ) THEN
450  RETURN
451  END IF
452 *
453 * Quick return if possible
454 *
455  IF( n.EQ.0 )
456  $ RETURN
457 *
458 * Get machine constants
459 *
460  eps = dlamch( 'P' )
461  smlnum = dlamch( 'S' )
462  bignum = one / smlnum
463  CALL dlabad( smlnum, bignum )
464  smlnum = sqrt( smlnum ) / eps
465  bignum = one / smlnum
466 *
467 * Scale A if max element outside range [SMLNUM,BIGNUM]
468 *
469  icond = 0
470  anrm = zlange( 'M', n, n, a, lda, dum )
471  scalea = .false.
472  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
473  scalea = .true.
474  cscale = smlnum
475  ELSE IF( anrm.GT.bignum ) THEN
476  scalea = .true.
477  cscale = bignum
478  END IF
479  IF( scalea )
480  $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
481 *
482 * Balance the matrix and compute ABNRM
483 *
484  CALL zgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
485  abnrm = zlange( '1', n, n, a, lda, dum )
486  IF( scalea ) THEN
487  dum( 1 ) = abnrm
488  CALL dlascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
489  abnrm = dum( 1 )
490  END IF
491 *
492 * Reduce to upper Hessenberg form
493 * (CWorkspace: need 2*N, prefer N+N*NB)
494 * (RWorkspace: none)
495 *
496  itau = 1
497  iwrk = itau + n
498  CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
499  $ lwork-iwrk+1, ierr )
500 *
501  IF( wantvl ) THEN
502 *
503 * Want left eigenvectors
504 * Copy Householder vectors to VL
505 *
506  side = 'L'
507  CALL zlacpy( 'L', n, n, a, lda, vl, ldvl )
508 *
509 * Generate unitary matrix in VL
510 * (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
511 * (RWorkspace: none)
512 *
513  CALL zunghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
514  $ lwork-iwrk+1, ierr )
515 *
516 * Perform QR iteration, accumulating Schur vectors in VL
517 * (CWorkspace: need 1, prefer HSWORK (see comments) )
518 * (RWorkspace: none)
519 *
520  iwrk = itau
521  CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vl, ldvl,
522  $ work( iwrk ), lwork-iwrk+1, info )
523 *
524  IF( wantvr ) THEN
525 *
526 * Want left and right eigenvectors
527 * Copy Schur vectors to VR
528 *
529  side = 'B'
530  CALL zlacpy( 'F', n, n, vl, ldvl, vr, ldvr )
531  END IF
532 *
533  ELSE IF( wantvr ) THEN
534 *
535 * Want right eigenvectors
536 * Copy Householder vectors to VR
537 *
538  side = 'R'
539  CALL zlacpy( 'L', n, n, a, lda, vr, ldvr )
540 *
541 * Generate unitary matrix in VR
542 * (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
543 * (RWorkspace: none)
544 *
545  CALL zunghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
546  $ lwork-iwrk+1, ierr )
547 *
548 * Perform QR iteration, accumulating Schur vectors in VR
549 * (CWorkspace: need 1, prefer HSWORK (see comments) )
550 * (RWorkspace: none)
551 *
552  iwrk = itau
553  CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vr, ldvr,
554  $ work( iwrk ), lwork-iwrk+1, info )
555 *
556  ELSE
557 *
558 * Compute eigenvalues only
559 * If condition numbers desired, compute Schur form
560 *
561  IF( wntsnn ) THEN
562  job = 'E'
563  ELSE
564  job = 'S'
565  END IF
566 *
567 * (CWorkspace: need 1, prefer HSWORK (see comments) )
568 * (RWorkspace: none)
569 *
570  iwrk = itau
571  CALL zhseqr( job, 'N', n, ilo, ihi, a, lda, w, vr, ldvr,
572  $ work( iwrk ), lwork-iwrk+1, info )
573  END IF
574 *
575 * If INFO .NE. 0 from ZHSEQR, then quit
576 *
577  IF( info.NE.0 )
578  $ GO TO 50
579 *
580  IF( wantvl .OR. wantvr ) THEN
581 *
582 * Compute left and/or right eigenvectors
583 * (CWorkspace: need 2*N, prefer N + 2*N*NB)
584 * (RWorkspace: need N)
585 *
586  CALL ztrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
587  $ n, nout, work( iwrk ), lwork-iwrk+1,
588  $ rwork, n, ierr )
589  END IF
590 *
591 * Compute condition numbers if desired
592 * (CWorkspace: need N*N+2*N unless SENSE = 'E')
593 * (RWorkspace: need 2*N unless SENSE = 'E')
594 *
595  IF( .NOT.wntsnn ) THEN
596  CALL ztrsna( sense, 'A', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
597  $ rconde, rcondv, n, nout, work( iwrk ), n, rwork,
598  $ icond )
599  END IF
600 *
601  IF( wantvl ) THEN
602 *
603 * Undo balancing of left eigenvectors
604 *
605  CALL zgebak( balanc, 'L', n, ilo, ihi, scale, n, vl, ldvl,
606  $ ierr )
607 *
608 * Normalize left eigenvectors and make largest component real
609 *
610  DO 20 i = 1, n
611  scl = one / dznrm2( n, vl( 1, i ), 1 )
612  CALL zdscal( n, scl, vl( 1, i ), 1 )
613  DO 10 k = 1, n
614  rwork( k ) = dble( vl( k, i ) )**2 +
615  $ aimag( vl( k, i ) )**2
616  10 CONTINUE
617  k = idamax( n, rwork, 1 )
618  tmp = conjg( vl( k, i ) ) / sqrt( rwork( k ) )
619  CALL zscal( n, tmp, vl( 1, i ), 1 )
620  vl( k, i ) = dcmplx( dble( vl( k, i ) ), zero )
621  20 CONTINUE
622  END IF
623 *
624  IF( wantvr ) THEN
625 *
626 * Undo balancing of right eigenvectors
627 *
628  CALL zgebak( balanc, 'R', n, ilo, ihi, scale, n, vr, ldvr,
629  $ ierr )
630 *
631 * Normalize right eigenvectors and make largest component real
632 *
633  DO 40 i = 1, n
634  scl = one / dznrm2( n, vr( 1, i ), 1 )
635  CALL zdscal( n, scl, vr( 1, i ), 1 )
636  DO 30 k = 1, n
637  rwork( k ) = dble( vr( k, i ) )**2 +
638  $ aimag( vr( k, i ) )**2
639  30 CONTINUE
640  k = idamax( n, rwork, 1 )
641  tmp = conjg( vr( k, i ) ) / sqrt( rwork( k ) )
642  CALL zscal( n, tmp, vr( 1, i ), 1 )
643  vr( k, i ) = dcmplx( dble( vr( k, i ) ), zero )
644  40 CONTINUE
645  END IF
646 *
647 * Undo scaling if necessary
648 *
649  50 CONTINUE
650  IF( scalea ) THEN
651  CALL zlascl( 'G', 0, 0, cscale, anrm, n-info, 1, w( info+1 ),
652  $ max( n-info, 1 ), ierr )
653  IF( info.EQ.0 ) THEN
654  IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
655  $ CALL dlascl( 'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
656  $ ierr )
657  ELSE
658  CALL zlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, w, n, ierr )
659  END IF
660  END IF
661 *
662  work( 1 ) = maxwrk
663  RETURN
664 *
665 * End of ZGEEVX
666 *
667  END
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zgebak(JOB, SIDE, N, ILO, IHI, SCALE, M, V, LDV, INFO)
ZGEBAK
Definition: zgebak.f:133
subroutine zgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZGEHRD
Definition: zgehrd.f:169
subroutine dlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: dlascl.f:145
subroutine zunghr(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGHR
Definition: zunghr.f:128
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zgebal(JOB, N, A, LDA, ILO, IHI, SCALE, INFO)
ZGEBAL
Definition: zgebal.f:162
subroutine zhseqr(JOB, COMPZ, N, ILO, IHI, H, LDH, W, Z, LDZ, WORK, LWORK, INFO)
ZHSEQR
Definition: zhseqr.f:301
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
subroutine zgeevx(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, W, VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE, RCONDV, WORK, LWORK, RWORK, INFO)
ZGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices ...
Definition: zgeevx.f:289
subroutine zlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: zlascl.f:145
subroutine ztrsna(JOB, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, S, SEP, MM, M, WORK, LDWORK, RWORK, INFO)
ZTRSNA
Definition: ztrsna.f:251
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine ztrevc3(SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, LWORK, RWORK, LRWORK, INFO)
ZTREVC3
Definition: ztrevc3.f:248