LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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*
280* @precisions fortran z -> c
281*
282*> \ingroup geevx
283*
284* =====================================================================
285 SUBROUTINE zgeevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, W, VL,
286 $ LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE,
287 $ RCONDV, WORK, LWORK, RWORK, INFO )
288 implicit none
289*
290* -- LAPACK driver routine --
291* -- LAPACK is a software package provided by Univ. of Tennessee, --
292* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
293*
294* .. Scalar Arguments ..
295 CHARACTER BALANC, JOBVL, JOBVR, SENSE
296 INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
297 DOUBLE PRECISION ABNRM
298* ..
299* .. Array Arguments ..
300 DOUBLE PRECISION RCONDE( * ), RCONDV( * ), RWORK( * ),
301 $ SCALE( * )
302 COMPLEX*16 A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ),
303 $ W( * ), WORK( * )
304* ..
305*
306* =====================================================================
307*
308* .. Parameters ..
309 DOUBLE PRECISION ZERO, ONE
310 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
311* ..
312* .. Local Scalars ..
313 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
314 $ WNTSNN, WNTSNV
315 CHARACTER JOB, SIDE
316 INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
317 $ lwork_trevc, maxwrk, minwrk, nout
318 DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, SCL, SMLNUM
319 COMPLEX*16 TMP
320* ..
321* .. Local Arrays ..
322 LOGICAL SELECT( 1 )
323 DOUBLE PRECISION DUM( 1 )
324* ..
325* .. External Subroutines ..
326 EXTERNAL dlascl, xerbla, zdscal, zgebak, zgebal, zgehrd,
328 $ zunghr
329* ..
330* .. External Functions ..
331 LOGICAL LSAME
332 INTEGER IDAMAX, ILAENV
333 DOUBLE PRECISION DLAMCH, DZNRM2, ZLANGE
334 EXTERNAL lsame, idamax, ilaenv, dlamch, dznrm2, zlange
335* ..
336* .. Intrinsic Functions ..
337 INTRINSIC dble, dcmplx, conjg, aimag, max, sqrt
338* ..
339* .. Executable Statements ..
340*
341* Test the input arguments
342*
343 info = 0
344 lquery = ( lwork.EQ.-1 )
345 wantvl = lsame( jobvl, 'V' )
346 wantvr = lsame( jobvr, 'V' )
347 wntsnn = lsame( sense, 'N' )
348 wntsne = lsame( sense, 'E' )
349 wntsnv = lsame( sense, 'V' )
350 wntsnb = lsame( sense, 'B' )
351 IF( .NOT.( lsame( balanc, 'N' ) .OR. lsame( balanc, 'S' ) .OR.
352 $ lsame( balanc, 'P' ) .OR. lsame( balanc, 'B' ) ) ) THEN
353 info = -1
354 ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
355 info = -2
356 ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
357 info = -3
358 ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
359 $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
360 $ wantvr ) ) ) THEN
361 info = -4
362 ELSE IF( n.LT.0 ) THEN
363 info = -5
364 ELSE IF( lda.LT.max( 1, n ) ) THEN
365 info = -7
366 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
367 info = -10
368 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
369 info = -12
370 END IF
371*
372* Compute workspace
373* (Note: Comments in the code beginning "Workspace:" describe the
374* minimal amount of workspace needed at that point in the code,
375* as well as the preferred amount for good performance.
376* CWorkspace refers to complex workspace, and RWorkspace to real
377* workspace. NB refers to the optimal block size for the
378* immediately following subroutine, as returned by ILAENV.
379* HSWORK refers to the workspace preferred by ZHSEQR, as
380* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
381* the worst case.)
382*
383 IF( info.EQ.0 ) THEN
384 IF( n.EQ.0 ) THEN
385 minwrk = 1
386 maxwrk = 1
387 ELSE
388 maxwrk = n + n*ilaenv( 1, 'ZGEHRD', ' ', n, 1, n, 0 )
389*
390 IF( wantvl ) THEN
391 CALL ztrevc3( 'L', 'B', SELECT, n, a, lda,
392 $ vl, ldvl, vr, ldvr,
393 $ n, nout, work, -1, rwork, -1, ierr )
394 lwork_trevc = int( work(1) )
395 maxwrk = max( maxwrk, lwork_trevc )
396 CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vl, ldvl,
397 $ work, -1, info )
398 ELSE IF( wantvr ) THEN
399 CALL ztrevc3( 'R', 'B', SELECT, n, a, lda,
400 $ vl, ldvl, vr, ldvr,
401 $ n, nout, work, -1, rwork, -1, ierr )
402 lwork_trevc = int( work(1) )
403 maxwrk = max( maxwrk, lwork_trevc )
404 CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vr, ldvr,
405 $ work, -1, info )
406 ELSE
407 IF( wntsnn ) THEN
408 CALL zhseqr( 'E', 'N', n, 1, n, a, lda, w, vr, ldvr,
409 $ work, -1, info )
410 ELSE
411 CALL zhseqr( 'S', 'N', n, 1, n, a, lda, w, vr, ldvr,
412 $ work, -1, info )
413 END IF
414 END IF
415 hswork = int( work(1) )
416*
417 IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) ) THEN
418 minwrk = 2*n
419 IF( .NOT.( wntsnn .OR. wntsne ) )
420 $ minwrk = max( minwrk, n*n + 2*n )
421 maxwrk = max( maxwrk, hswork )
422 IF( .NOT.( wntsnn .OR. wntsne ) )
423 $ maxwrk = max( maxwrk, n*n + 2*n )
424 ELSE
425 minwrk = 2*n
426 IF( .NOT.( wntsnn .OR. wntsne ) )
427 $ minwrk = max( minwrk, n*n + 2*n )
428 maxwrk = max( maxwrk, hswork )
429 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1, 'ZUNGHR',
430 $ ' ', n, 1, n, -1 ) )
431 IF( .NOT.( wntsnn .OR. wntsne ) )
432 $ maxwrk = max( maxwrk, n*n + 2*n )
433 maxwrk = max( maxwrk, 2*n )
434 END IF
435 maxwrk = max( maxwrk, minwrk )
436 END IF
437 work( 1 ) = maxwrk
438*
439 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
440 info = -20
441 END IF
442 END IF
443*
444 IF( info.NE.0 ) THEN
445 CALL xerbla( 'ZGEEVX', -info )
446 RETURN
447 ELSE IF( lquery ) THEN
448 RETURN
449 END IF
450*
451* Quick return if possible
452*
453 IF( n.EQ.0 )
454 $ RETURN
455*
456* Get machine constants
457*
458 eps = dlamch( 'P' )
459 smlnum = dlamch( 'S' )
460 bignum = one / smlnum
461 smlnum = sqrt( smlnum ) / eps
462 bignum = one / smlnum
463*
464* Scale A if max element outside range [SMLNUM,BIGNUM]
465*
466 icond = 0
467 anrm = zlange( 'M', n, n, a, lda, dum )
468 scalea = .false.
469 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
470 scalea = .true.
471 cscale = smlnum
472 ELSE IF( anrm.GT.bignum ) THEN
473 scalea = .true.
474 cscale = bignum
475 END IF
476 IF( scalea )
477 $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
478*
479* Balance the matrix and compute ABNRM
480*
481 CALL zgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
482 abnrm = zlange( '1', n, n, a, lda, dum )
483 IF( scalea ) THEN
484 dum( 1 ) = abnrm
485 CALL dlascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
486 abnrm = dum( 1 )
487 END IF
488*
489* Reduce to upper Hessenberg form
490* (CWorkspace: need 2*N, prefer N+N*NB)
491* (RWorkspace: none)
492*
493 itau = 1
494 iwrk = itau + n
495 CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
496 $ lwork-iwrk+1, ierr )
497*
498 IF( wantvl ) THEN
499*
500* Want left eigenvectors
501* Copy Householder vectors to VL
502*
503 side = 'L'
504 CALL zlacpy( 'L', n, n, a, lda, vl, ldvl )
505*
506* Generate unitary matrix in VL
507* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
508* (RWorkspace: none)
509*
510 CALL zunghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
511 $ lwork-iwrk+1, ierr )
512*
513* Perform QR iteration, accumulating Schur vectors in VL
514* (CWorkspace: need 1, prefer HSWORK (see comments) )
515* (RWorkspace: none)
516*
517 iwrk = itau
518 CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vl, ldvl,
519 $ work( iwrk ), lwork-iwrk+1, info )
520*
521 IF( wantvr ) THEN
522*
523* Want left and right eigenvectors
524* Copy Schur vectors to VR
525*
526 side = 'B'
527 CALL zlacpy( 'F', n, n, vl, ldvl, vr, ldvr )
528 END IF
529*
530 ELSE IF( wantvr ) THEN
531*
532* Want right eigenvectors
533* Copy Householder vectors to VR
534*
535 side = 'R'
536 CALL zlacpy( 'L', n, n, a, lda, vr, ldvr )
537*
538* Generate unitary matrix in VR
539* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
540* (RWorkspace: none)
541*
542 CALL zunghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
543 $ lwork-iwrk+1, ierr )
544*
545* Perform QR iteration, accumulating Schur vectors in VR
546* (CWorkspace: need 1, prefer HSWORK (see comments) )
547* (RWorkspace: none)
548*
549 iwrk = itau
550 CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vr, ldvr,
551 $ work( iwrk ), lwork-iwrk+1, info )
552*
553 ELSE
554*
555* Compute eigenvalues only
556* If condition numbers desired, compute Schur form
557*
558 IF( wntsnn ) THEN
559 job = 'E'
560 ELSE
561 job = 'S'
562 END IF
563*
564* (CWorkspace: need 1, prefer HSWORK (see comments) )
565* (RWorkspace: none)
566*
567 iwrk = itau
568 CALL zhseqr( job, 'N', n, ilo, ihi, a, lda, w, vr, ldvr,
569 $ work( iwrk ), lwork-iwrk+1, info )
570 END IF
571*
572* If INFO .NE. 0 from ZHSEQR, then quit
573*
574 IF( info.NE.0 )
575 $ GO TO 50
576*
577 IF( wantvl .OR. wantvr ) THEN
578*
579* Compute left and/or right eigenvectors
580* (CWorkspace: need 2*N, prefer N + 2*N*NB)
581* (RWorkspace: need N)
582*
583 CALL ztrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
584 $ n, nout, work( iwrk ), lwork-iwrk+1,
585 $ rwork, n, ierr )
586 END IF
587*
588* Compute condition numbers if desired
589* (CWorkspace: need N*N+2*N unless SENSE = 'E')
590* (RWorkspace: need 2*N unless SENSE = 'E')
591*
592 IF( .NOT.wntsnn ) THEN
593 CALL ztrsna( sense, 'A', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
594 $ rconde, rcondv, n, nout, work( iwrk ), n, rwork,
595 $ icond )
596 END IF
597*
598 IF( wantvl ) THEN
599*
600* Undo balancing of left eigenvectors
601*
602 CALL zgebak( balanc, 'L', n, ilo, ihi, scale, n, vl, ldvl,
603 $ ierr )
604*
605* Normalize left eigenvectors and make largest component real
606*
607 DO 20 i = 1, n
608 scl = one / dznrm2( n, vl( 1, i ), 1 )
609 CALL zdscal( n, scl, vl( 1, i ), 1 )
610 DO 10 k = 1, n
611 rwork( k ) = dble( vl( k, i ) )**2 +
612 $ aimag( vl( k, i ) )**2
613 10 CONTINUE
614 k = idamax( n, rwork, 1 )
615 tmp = conjg( vl( k, i ) ) / sqrt( rwork( k ) )
616 CALL zscal( n, tmp, vl( 1, i ), 1 )
617 vl( k, i ) = dcmplx( dble( vl( k, i ) ), zero )
618 20 CONTINUE
619 END IF
620*
621 IF( wantvr ) THEN
622*
623* Undo balancing of right eigenvectors
624*
625 CALL zgebak( balanc, 'R', n, ilo, ihi, scale, n, vr, ldvr,
626 $ ierr )
627*
628* Normalize right eigenvectors and make largest component real
629*
630 DO 40 i = 1, n
631 scl = one / dznrm2( n, vr( 1, i ), 1 )
632 CALL zdscal( n, scl, vr( 1, i ), 1 )
633 DO 30 k = 1, n
634 rwork( k ) = dble( vr( k, i ) )**2 +
635 $ aimag( vr( k, i ) )**2
636 30 CONTINUE
637 k = idamax( n, rwork, 1 )
638 tmp = conjg( vr( k, i ) ) / sqrt( rwork( k ) )
639 CALL zscal( n, tmp, vr( 1, i ), 1 )
640 vr( k, i ) = dcmplx( dble( vr( k, i ) ), zero )
641 40 CONTINUE
642 END IF
643*
644* Undo scaling if necessary
645*
646 50 CONTINUE
647 IF( scalea ) THEN
648 CALL zlascl( 'G', 0, 0, cscale, anrm, n-info, 1, w( info+1 ),
649 $ max( n-info, 1 ), ierr )
650 IF( info.EQ.0 ) THEN
651 IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
652 $ CALL dlascl( 'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
653 $ ierr )
654 ELSE
655 CALL zlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, w, n, ierr )
656 END IF
657 END IF
658*
659 work( 1 ) = maxwrk
660 RETURN
661*
662* End of ZGEEVX
663*
664 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
ZGEBAK
Definition zgebak.f:131
subroutine zgebal(job, n, a, lda, ilo, ihi, scale, info)
ZGEBAL
Definition zgebal.f:165
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:288
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
Definition zgehrd.f:167
subroutine zhseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
ZHSEQR
Definition zhseqr.f:299
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
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:143
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:143
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zscal(n, za, zx, incx)
ZSCAL
Definition zscal.f:78
subroutine ztrevc3(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, lwork, rwork, lrwork, info)
ZTREVC3
Definition ztrevc3.f:244
subroutine ztrsna(job, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, s, sep, mm, m, work, ldwork, rwork, info)
ZTRSNA
Definition ztrsna.f:249
subroutine zunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZUNGHR
Definition zunghr.f:126