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