LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cggevx.f
Go to the documentation of this file.
1*> \brief <b> CGGEVX 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 CGGEVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cggevx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cggevx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cggevx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
20* ALPHA, BETA, VL, LDVL, VR, LDVR, ILO, IHI,
21* LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV,
22* WORK, LWORK, RWORK, IWORK, BWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER BALANC, JOBVL, JOBVR, SENSE
26* INTEGER IHI, ILO, INFO, LDA, LDB, LDVL, LDVR, LWORK, N
27* REAL ABNRM, BBNRM
28* ..
29* .. Array Arguments ..
30* LOGICAL BWORK( * )
31* INTEGER IWORK( * )
32* REAL LSCALE( * ), RCONDE( * ), RCONDV( * ),
33* $ RSCALE( * ), RWORK( * )
34* COMPLEX A( LDA, * ), ALPHA( * ), B( LDB, * ),
35* $ BETA( * ), VL( LDVL, * ), VR( LDVR, * ),
36* $ WORK( * )
37* ..
38*
39*
40*> \par Purpose:
41* =============
42*>
43*> \verbatim
44*>
45*> CGGEVX computes for a pair of N-by-N complex nonsymmetric matrices
46*> (A,B) the generalized eigenvalues, and optionally, the left and/or
47*> right generalized eigenvectors.
48*>
49*> Optionally, it also computes a balancing transformation to improve
50*> the conditioning of the eigenvalues and eigenvectors (ILO, IHI,
51*> LSCALE, RSCALE, ABNRM, and BBNRM), reciprocal condition numbers for
52*> the eigenvalues (RCONDE), and reciprocal condition numbers for the
53*> right eigenvectors (RCONDV).
54*>
55*> A generalized eigenvalue for a pair of matrices (A,B) is a scalar
56*> lambda or a ratio alpha/beta = lambda, such that A - lambda*B is
57*> singular. It is usually represented as the pair (alpha,beta), as
58*> there is a reasonable interpretation for beta=0, and even for both
59*> being zero.
60*>
61*> The right eigenvector v(j) corresponding to the eigenvalue lambda(j)
62*> of (A,B) satisfies
63*> A * v(j) = lambda(j) * B * v(j) .
64*> The left eigenvector u(j) corresponding to the eigenvalue lambda(j)
65*> of (A,B) satisfies
66*> u(j)**H * A = lambda(j) * u(j)**H * B.
67*> where u(j)**H is the conjugate-transpose of u(j).
68*>
69*> \endverbatim
70*
71* Arguments:
72* ==========
73*
74*> \param[in] BALANC
75*> \verbatim
76*> BALANC is CHARACTER*1
77*> Specifies the balance option to be performed:
78*> = 'N': do not diagonally scale or permute;
79*> = 'P': permute only;
80*> = 'S': scale only;
81*> = 'B': both permute and scale.
82*> Computed reciprocal condition numbers will be for the
83*> matrices after permuting and/or balancing. Permuting does
84*> not change condition numbers (in exact arithmetic), but
85*> balancing does.
86*> \endverbatim
87*>
88*> \param[in] JOBVL
89*> \verbatim
90*> JOBVL is CHARACTER*1
91*> = 'N': do not compute the left generalized eigenvectors;
92*> = 'V': compute the left generalized eigenvectors.
93*> \endverbatim
94*>
95*> \param[in] JOBVR
96*> \verbatim
97*> JOBVR is CHARACTER*1
98*> = 'N': do not compute the right generalized eigenvectors;
99*> = 'V': compute the right generalized eigenvectors.
100*> \endverbatim
101*>
102*> \param[in] SENSE
103*> \verbatim
104*> SENSE is CHARACTER*1
105*> Determines which reciprocal condition numbers are computed.
106*> = 'N': none are computed;
107*> = 'E': computed for eigenvalues only;
108*> = 'V': computed for eigenvectors only;
109*> = 'B': computed for eigenvalues and eigenvectors.
110*> \endverbatim
111*>
112*> \param[in] N
113*> \verbatim
114*> N is INTEGER
115*> The order of the matrices A, B, VL, and VR. N >= 0.
116*> \endverbatim
117*>
118*> \param[in,out] A
119*> \verbatim
120*> A is COMPLEX array, dimension (LDA, N)
121*> On entry, the matrix A in the pair (A,B).
122*> On exit, A has been overwritten. If JOBVL='V' or JOBVR='V'
123*> or both, then A contains the first part of the complex Schur
124*> form of the "balanced" versions of the input A and B.
125*> \endverbatim
126*>
127*> \param[in] LDA
128*> \verbatim
129*> LDA is INTEGER
130*> The leading dimension of A. LDA >= max(1,N).
131*> \endverbatim
132*>
133*> \param[in,out] B
134*> \verbatim
135*> B is COMPLEX array, dimension (LDB, N)
136*> On entry, the matrix B in the pair (A,B).
137*> On exit, B has been overwritten. If JOBVL='V' or JOBVR='V'
138*> or both, then B contains the second part of the complex
139*> Schur form of the "balanced" versions of the input A and B.
140*> \endverbatim
141*>
142*> \param[in] LDB
143*> \verbatim
144*> LDB is INTEGER
145*> The leading dimension of B. LDB >= max(1,N).
146*> \endverbatim
147*>
148*> \param[out] ALPHA
149*> \verbatim
150*> ALPHA is COMPLEX array, dimension (N)
151*> \endverbatim
152*>
153*> \param[out] BETA
154*> \verbatim
155*> BETA is COMPLEX array, dimension (N)
156*> On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the generalized
157*> eigenvalues.
158*>
159*> Note: the quotient ALPHA(j)/BETA(j) ) may easily over- or
160*> underflow, and BETA(j) may even be zero. Thus, the user
161*> should avoid naively computing the ratio ALPHA/BETA.
162*> However, ALPHA will be always less than and usually
163*> comparable with norm(A) in magnitude, and BETA always less
164*> than and usually comparable with norm(B).
165*> \endverbatim
166*>
167*> \param[out] VL
168*> \verbatim
169*> VL is COMPLEX array, dimension (LDVL,N)
170*> If JOBVL = 'V', the left generalized eigenvectors u(j) are
171*> stored one after another in the columns of VL, in the same
172*> order as their eigenvalues.
173*> Each eigenvector will be scaled so the largest component
174*> will have abs(real part) + abs(imag. part) = 1.
175*> Not referenced if JOBVL = 'N'.
176*> \endverbatim
177*>
178*> \param[in] LDVL
179*> \verbatim
180*> LDVL is INTEGER
181*> The leading dimension of the matrix VL. LDVL >= 1, and
182*> if JOBVL = 'V', LDVL >= N.
183*> \endverbatim
184*>
185*> \param[out] VR
186*> \verbatim
187*> VR is COMPLEX array, dimension (LDVR,N)
188*> If JOBVR = 'V', the right generalized eigenvectors v(j) are
189*> stored one after another in the columns of VR, in the same
190*> order as their eigenvalues.
191*> Each eigenvector will be scaled so the largest component
192*> will have abs(real part) + abs(imag. part) = 1.
193*> Not referenced if JOBVR = 'N'.
194*> \endverbatim
195*>
196*> \param[in] LDVR
197*> \verbatim
198*> LDVR is INTEGER
199*> The leading dimension of the matrix VR. LDVR >= 1, and
200*> if JOBVR = 'V', LDVR >= N.
201*> \endverbatim
202*>
203*> \param[out] ILO
204*> \verbatim
205*> ILO is INTEGER
206*> \endverbatim
207*>
208*> \param[out] IHI
209*> \verbatim
210*> IHI is INTEGER
211*> ILO and IHI are integer values such that on exit
212*> A(i,j) = 0 and B(i,j) = 0 if i > j and
213*> j = 1,...,ILO-1 or i = IHI+1,...,N.
214*> If BALANC = 'N' or 'S', ILO = 1 and IHI = N.
215*> \endverbatim
216*>
217*> \param[out] LSCALE
218*> \verbatim
219*> LSCALE is REAL array, dimension (N)
220*> Details of the permutations and scaling factors applied
221*> to the left side of A and B. If PL(j) is the index of the
222*> row interchanged with row j, and DL(j) is the scaling
223*> factor applied to row j, then
224*> LSCALE(j) = PL(j) for j = 1,...,ILO-1
225*> = DL(j) for j = ILO,...,IHI
226*> = PL(j) for j = IHI+1,...,N.
227*> The order in which the interchanges are made is N to IHI+1,
228*> then 1 to ILO-1.
229*> \endverbatim
230*>
231*> \param[out] RSCALE
232*> \verbatim
233*> RSCALE is REAL array, dimension (N)
234*> Details of the permutations and scaling factors applied
235*> to the right side of A and B. If PR(j) is the index of the
236*> column interchanged with column j, and DR(j) is the scaling
237*> factor applied to column j, then
238*> RSCALE(j) = PR(j) for j = 1,...,ILO-1
239*> = DR(j) for j = ILO,...,IHI
240*> = PR(j) for j = IHI+1,...,N
241*> The order in which the interchanges are made is N to IHI+1,
242*> then 1 to ILO-1.
243*> \endverbatim
244*>
245*> \param[out] ABNRM
246*> \verbatim
247*> ABNRM is REAL
248*> The one-norm of the balanced matrix A.
249*> \endverbatim
250*>
251*> \param[out] BBNRM
252*> \verbatim
253*> BBNRM is REAL
254*> The one-norm of the balanced matrix B.
255*> \endverbatim
256*>
257*> \param[out] RCONDE
258*> \verbatim
259*> RCONDE is REAL array, dimension (N)
260*> If SENSE = 'E' or 'B', the reciprocal condition numbers of
261*> the eigenvalues, stored in consecutive elements of the array.
262*> If SENSE = 'N' or 'V', RCONDE is not referenced.
263*> \endverbatim
264*>
265*> \param[out] RCONDV
266*> \verbatim
267*> RCONDV is REAL array, dimension (N)
268*> If SENSE = 'V' or 'B', the estimated reciprocal condition
269*> numbers of the eigenvectors, stored in consecutive elements
270*> of the array. If the eigenvalues cannot be reordered to
271*> compute RCONDV(j), RCONDV(j) is set to 0; this can only occur
272*> when the true value would be very small anyway.
273*> If SENSE = 'N' or 'E', RCONDV is not referenced.
274*> \endverbatim
275*>
276*> \param[out] WORK
277*> \verbatim
278*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
279*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
280*> \endverbatim
281*>
282*> \param[in] LWORK
283*> \verbatim
284*> LWORK is INTEGER
285*> The dimension of the array WORK. LWORK >= max(1,2*N).
286*> If SENSE = 'E', LWORK >= max(1,4*N).
287*> If SENSE = 'V' or 'B', LWORK >= max(1,2*N*N+2*N).
288*>
289*> If LWORK = -1, then a workspace query is assumed; the routine
290*> only calculates the optimal size of the WORK array, returns
291*> this value as the first entry of the WORK array, and no error
292*> message related to LWORK is issued by XERBLA.
293*> \endverbatim
294*>
295*> \param[out] RWORK
296*> \verbatim
297*> RWORK is REAL array, dimension (lrwork)
298*> lrwork must be at least max(1,6*N) if BALANC = 'S' or 'B',
299*> and at least max(1,2*N) otherwise.
300*> Real workspace.
301*> \endverbatim
302*>
303*> \param[out] IWORK
304*> \verbatim
305*> IWORK is INTEGER array, dimension (N+2)
306*> If SENSE = 'E', IWORK is not referenced.
307*> \endverbatim
308*>
309*> \param[out] BWORK
310*> \verbatim
311*> BWORK is LOGICAL array, dimension (N)
312*> If SENSE = 'N', BWORK is not referenced.
313*> \endverbatim
314*>
315*> \param[out] INFO
316*> \verbatim
317*> INFO is INTEGER
318*> = 0: successful exit
319*> < 0: if INFO = -i, the i-th argument had an illegal value.
320*> = 1,...,N:
321*> The QZ iteration failed. No eigenvectors have been
322*> calculated, but ALPHA(j) and BETA(j) should be correct
323*> for j=INFO+1,...,N.
324*> > N: =N+1: other than QZ iteration failed in CHGEQZ.
325*> =N+2: error return from CTGEVC.
326*> \endverbatim
327*
328* Authors:
329* ========
330*
331*> \author Univ. of Tennessee
332*> \author Univ. of California Berkeley
333*> \author Univ. of Colorado Denver
334*> \author NAG Ltd.
335*
336*> \ingroup ggevx
337*
338*> \par Further Details:
339* =====================
340*>
341*> \verbatim
342*>
343*> Balancing a matrix pair (A,B) includes, first, permuting rows and
344*> columns to isolate eigenvalues, second, applying diagonal similarity
345*> transformation to the rows and columns to make the rows and columns
346*> as close in norm as possible. The computed reciprocal condition
347*> numbers correspond to the balanced matrix. Permuting rows and columns
348*> will not change the condition numbers (in exact arithmetic) but
349*> diagonal scaling will. For further explanation of balancing, see
350*> section 4.11.1.2 of LAPACK Users' Guide.
351*>
352*> An approximate error bound on the chordal distance between the i-th
353*> computed generalized eigenvalue w and the corresponding exact
354*> eigenvalue lambda is
355*>
356*> chord(w, lambda) <= EPS * norm(ABNRM, BBNRM) / RCONDE(I)
357*>
358*> An approximate error bound for the angle between the i-th computed
359*> eigenvector VL(i) or VR(i) is given by
360*>
361*> EPS * norm(ABNRM, BBNRM) / DIF(i).
362*>
363*> For further explanation of the reciprocal condition numbers RCONDE
364*> and RCONDV, see section 4.11 of LAPACK User's Guide.
365*> \endverbatim
366*>
367* =====================================================================
368 SUBROUTINE cggevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B,
369 $ LDB,
370 $ ALPHA, BETA, VL, LDVL, VR, LDVR, ILO, IHI,
371 $ LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV,
372 $ WORK, LWORK, RWORK, IWORK, BWORK, INFO )
373*
374* -- LAPACK driver routine --
375* -- LAPACK is a software package provided by Univ. of Tennessee, --
376* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
377*
378* .. Scalar Arguments ..
379 CHARACTER BALANC, JOBVL, JOBVR, SENSE
380 INTEGER IHI, ILO, INFO, LDA, LDB, LDVL, LDVR, LWORK, N
381 REAL ABNRM, BBNRM
382* ..
383* .. Array Arguments ..
384 LOGICAL BWORK( * )
385 INTEGER IWORK( * )
386 REAL LSCALE( * ), RCONDE( * ), RCONDV( * ),
387 $ RSCALE( * ), RWORK( * )
388 COMPLEX A( LDA, * ), ALPHA( * ), B( LDB, * ),
389 $ beta( * ), vl( ldvl, * ), vr( ldvr, * ),
390 $ work( * )
391* ..
392*
393* =====================================================================
394*
395* .. Parameters ..
396 REAL ZERO, ONE
397 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
398 COMPLEX CZERO, CONE
399 PARAMETER ( CZERO = ( 0.0e+0, 0.0e+0 ),
400 $ cone = ( 1.0e+0, 0.0e+0 ) )
401* ..
402* .. Local Scalars ..
403 LOGICAL ILASCL, ILBSCL, ILV, ILVL, ILVR, LQUERY, NOSCL,
404 $ WANTSB, WANTSE, WANTSN, WANTSV
405 CHARACTER CHTEMP
406 INTEGER I, ICOLS, IERR, IJOBVL, IJOBVR, IN, IROWS,
407 $ itau, iwrk, iwrk1, j, jc, jr, m, maxwrk, minwrk
408 REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
409 $ smlnum, temp
410 COMPLEX X
411* ..
412* .. Local Arrays ..
413 LOGICAL LDUMMA( 1 )
414* ..
415* .. External Subroutines ..
416 EXTERNAL cgeqrf, cggbak, cggbal, cgghrd, chgeqz,
417 $ clacpy,
419 $ slascl, xerbla
420* ..
421* .. External Functions ..
422 LOGICAL LSAME
423 INTEGER ILAENV
424 REAL CLANGE, SLAMCH, SROUNDUP_LWORK
425 EXTERNAL LSAME, ILAENV, CLANGE, SLAMCH,
426 $ sroundup_lwork
427* ..
428* .. Intrinsic Functions ..
429 INTRINSIC abs, aimag, max, real, sqrt
430* ..
431* .. Statement Functions ..
432 REAL ABS1
433* ..
434* .. Statement Function definitions ..
435 abs1( x ) = abs( real( x ) ) + abs( aimag( x ) )
436* ..
437* .. Executable Statements ..
438*
439* Decode the input arguments
440*
441 IF( lsame( jobvl, 'N' ) ) THEN
442 ijobvl = 1
443 ilvl = .false.
444 ELSE IF( lsame( jobvl, 'V' ) ) THEN
445 ijobvl = 2
446 ilvl = .true.
447 ELSE
448 ijobvl = -1
449 ilvl = .false.
450 END IF
451*
452 IF( lsame( jobvr, 'N' ) ) THEN
453 ijobvr = 1
454 ilvr = .false.
455 ELSE IF( lsame( jobvr, 'V' ) ) THEN
456 ijobvr = 2
457 ilvr = .true.
458 ELSE
459 ijobvr = -1
460 ilvr = .false.
461 END IF
462 ilv = ilvl .OR. ilvr
463*
464 noscl = lsame( balanc, 'N' ) .OR. lsame( balanc, 'P' )
465 wantsn = lsame( sense, 'N' )
466 wantse = lsame( sense, 'E' )
467 wantsv = lsame( sense, 'V' )
468 wantsb = lsame( sense, 'B' )
469*
470* Test the input arguments
471*
472 info = 0
473 lquery = ( lwork.EQ.-1 )
474 IF( .NOT.( noscl .OR. lsame( balanc,'S' ) .OR.
475 $ lsame( balanc, 'B' ) ) ) THEN
476 info = -1
477 ELSE IF( ijobvl.LE.0 ) THEN
478 info = -2
479 ELSE IF( ijobvr.LE.0 ) THEN
480 info = -3
481 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsb .OR. wantsv ) )
482 $ THEN
483 info = -4
484 ELSE IF( n.LT.0 ) THEN
485 info = -5
486 ELSE IF( lda.LT.max( 1, n ) ) THEN
487 info = -7
488 ELSE IF( ldb.LT.max( 1, n ) ) THEN
489 info = -9
490 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
491 info = -13
492 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
493 info = -15
494 END IF
495*
496* Compute workspace
497* (Note: Comments in the code beginning "Workspace:" describe the
498* minimal amount of workspace needed at that point in the code,
499* as well as the preferred amount for good performance.
500* NB refers to the optimal block size for the immediately
501* following subroutine, as returned by ILAENV. The workspace is
502* computed assuming ILO = 1 and IHI = N, the worst case.)
503*
504 IF( info.EQ.0 ) THEN
505 IF( n.EQ.0 ) THEN
506 minwrk = 1
507 maxwrk = 1
508 ELSE
509 minwrk = 2*n
510 IF( wantse ) THEN
511 minwrk = 4*n
512 ELSE IF( wantsv .OR. wantsb ) THEN
513 minwrk = 2*n*( n + 1)
514 END IF
515 maxwrk = minwrk
516 maxwrk = max( maxwrk,
517 $ n + n*ilaenv( 1, 'CGEQRF', ' ', n, 1, n,
518 $ 0 ) )
519 maxwrk = max( maxwrk,
520 $ n + n*ilaenv( 1, 'CUNMQR', ' ', n, 1, n,
521 $ 0 ) )
522 IF( ilvl ) THEN
523 maxwrk = max( maxwrk, n +
524 $ n*ilaenv( 1, 'CUNGQR', ' ', n, 1, n,
525 $ 0 ) )
526 END IF
527 END IF
528 work( 1 ) = sroundup_lwork(maxwrk)
529*
530 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
531 info = -25
532 END IF
533 END IF
534*
535 IF( info.NE.0 ) THEN
536 CALL xerbla( 'CGGEVX', -info )
537 RETURN
538 ELSE IF( lquery ) THEN
539 RETURN
540 END IF
541*
542* Quick return if possible
543*
544 IF( n.EQ.0 )
545 $ RETURN
546*
547* Get machine constants
548*
549 eps = slamch( 'P' )
550 smlnum = slamch( 'S' )
551 bignum = one / smlnum
552 smlnum = sqrt( smlnum ) / eps
553 bignum = one / smlnum
554*
555* Scale A if max element outside range [SMLNUM,BIGNUM]
556*
557 anrm = clange( 'M', n, n, a, lda, rwork )
558 ilascl = .false.
559 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
560 anrmto = smlnum
561 ilascl = .true.
562 ELSE IF( anrm.GT.bignum ) THEN
563 anrmto = bignum
564 ilascl = .true.
565 END IF
566 IF( ilascl )
567 $ CALL clascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
568*
569* Scale B if max element outside range [SMLNUM,BIGNUM]
570*
571 bnrm = clange( 'M', n, n, b, ldb, rwork )
572 ilbscl = .false.
573 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
574 bnrmto = smlnum
575 ilbscl = .true.
576 ELSE IF( bnrm.GT.bignum ) THEN
577 bnrmto = bignum
578 ilbscl = .true.
579 END IF
580 IF( ilbscl )
581 $ CALL clascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
582*
583* Permute and/or balance the matrix pair (A,B)
584* (Real Workspace: need 6*N if BALANC = 'S' or 'B', 1 otherwise)
585*
586 CALL cggbal( balanc, n, a, lda, b, ldb, ilo, ihi, lscale,
587 $ rscale,
588 $ rwork, ierr )
589*
590* Compute ABNRM and BBNRM
591*
592 abnrm = clange( '1', n, n, a, lda, rwork( 1 ) )
593 IF( ilascl ) THEN
594 rwork( 1 ) = abnrm
595 CALL slascl( 'G', 0, 0, anrmto, anrm, 1, 1, rwork( 1 ), 1,
596 $ ierr )
597 abnrm = rwork( 1 )
598 END IF
599*
600 bbnrm = clange( '1', n, n, b, ldb, rwork( 1 ) )
601 IF( ilbscl ) THEN
602 rwork( 1 ) = bbnrm
603 CALL slascl( 'G', 0, 0, bnrmto, bnrm, 1, 1, rwork( 1 ), 1,
604 $ ierr )
605 bbnrm = rwork( 1 )
606 END IF
607*
608* Reduce B to triangular form (QR decomposition of B)
609* (Complex Workspace: need N, prefer N*NB )
610*
611 irows = ihi + 1 - ilo
612 IF( ilv .OR. .NOT.wantsn ) THEN
613 icols = n + 1 - ilo
614 ELSE
615 icols = irows
616 END IF
617 itau = 1
618 iwrk = itau + irows
619 CALL cgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
620 $ work( iwrk ), lwork+1-iwrk, ierr )
621*
622* Apply the unitary transformation to A
623* (Complex Workspace: need N, prefer N*NB)
624*
625 CALL cunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
626 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
627 $ lwork+1-iwrk, ierr )
628*
629* Initialize VL and/or VR
630* (Workspace: need N, prefer N*NB)
631*
632 IF( ilvl ) THEN
633 CALL claset( 'Full', n, n, czero, cone, vl, ldvl )
634 IF( irows.GT.1 ) THEN
635 CALL clacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
636 $ vl( ilo+1, ilo ), ldvl )
637 END IF
638 CALL cungqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
639 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
640 END IF
641*
642 IF( ilvr )
643 $ CALL claset( 'Full', n, n, czero, cone, vr, ldvr )
644*
645* Reduce to generalized Hessenberg form
646* (Workspace: none needed)
647*
648 IF( ilv .OR. .NOT.wantsn ) THEN
649*
650* Eigenvectors requested -- work on whole matrix.
651*
652 CALL cgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
653 $ ldvl, vr, ldvr, ierr )
654 ELSE
655 CALL cgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
656 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
657 END IF
658*
659* Perform QZ algorithm (Compute eigenvalues, and optionally, the
660* Schur forms and Schur vectors)
661* (Complex Workspace: need N)
662* (Real Workspace: need N)
663*
664 iwrk = itau
665 IF( ilv .OR. .NOT.wantsn ) THEN
666 chtemp = 'S'
667 ELSE
668 chtemp = 'E'
669 END IF
670*
671 CALL chgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
672 $ alpha, beta, vl, ldvl, vr, ldvr, work( iwrk ),
673 $ lwork+1-iwrk, rwork, ierr )
674 IF( ierr.NE.0 ) THEN
675 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
676 info = ierr
677 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
678 info = ierr - n
679 ELSE
680 info = n + 1
681 END IF
682 GO TO 90
683 END IF
684*
685* Compute Eigenvectors and estimate condition numbers if desired
686* CTGEVC: (Complex Workspace: need 2*N )
687* (Real Workspace: need 2*N )
688* CTGSNA: (Complex Workspace: need 2*N*N if SENSE='V' or 'B')
689* (Integer Workspace: need N+2 )
690*
691 IF( ilv .OR. .NOT.wantsn ) THEN
692 IF( ilv ) THEN
693 IF( ilvl ) THEN
694 IF( ilvr ) THEN
695 chtemp = 'B'
696 ELSE
697 chtemp = 'L'
698 END IF
699 ELSE
700 chtemp = 'R'
701 END IF
702*
703 CALL ctgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl,
704 $ ldvl, vr, ldvr, n, in, work( iwrk ), rwork,
705 $ ierr )
706 IF( ierr.NE.0 ) THEN
707 info = n + 2
708 GO TO 90
709 END IF
710 END IF
711*
712 IF( .NOT.wantsn ) THEN
713*
714* compute eigenvectors (CTGEVC) and estimate condition
715* numbers (CTGSNA). Note that the definition of the condition
716* number is not invariant under transformation (u,v) to
717* (Q*u, Z*v), where (u,v) are eigenvectors of the generalized
718* Schur form (S,T), Q and Z are orthogonal matrices. In order
719* to avoid using extra 2*N*N workspace, we have to
720* re-calculate eigenvectors and estimate the condition numbers
721* one at a time.
722*
723 DO 20 i = 1, n
724*
725 DO 10 j = 1, n
726 bwork( j ) = .false.
727 10 CONTINUE
728 bwork( i ) = .true.
729*
730 iwrk = n + 1
731 iwrk1 = iwrk + n
732*
733 IF( wantse .OR. wantsb ) THEN
734 CALL ctgevc( 'B', 'S', bwork, n, a, lda, b, ldb,
735 $ work( 1 ), n, work( iwrk ), n, 1, m,
736 $ work( iwrk1 ), rwork, ierr )
737 IF( ierr.NE.0 ) THEN
738 info = n + 2
739 GO TO 90
740 END IF
741 END IF
742*
743 CALL ctgsna( sense, 'S', bwork, n, a, lda, b, ldb,
744 $ work( 1 ), n, work( iwrk ), n, rconde( i ),
745 $ rcondv( i ), 1, m, work( iwrk1 ),
746 $ lwork-iwrk1+1, iwork, ierr )
747*
748 20 CONTINUE
749 END IF
750 END IF
751*
752* Undo balancing on VL and VR and normalization
753* (Workspace: none needed)
754*
755 IF( ilvl ) THEN
756 CALL cggbak( balanc, 'L', n, ilo, ihi, lscale, rscale, n,
757 $ vl,
758 $ ldvl, ierr )
759*
760 DO 50 jc = 1, n
761 temp = zero
762 DO 30 jr = 1, n
763 temp = max( temp, abs1( vl( jr, jc ) ) )
764 30 CONTINUE
765 IF( temp.LT.smlnum )
766 $ GO TO 50
767 temp = one / temp
768 DO 40 jr = 1, n
769 vl( jr, jc ) = vl( jr, jc )*temp
770 40 CONTINUE
771 50 CONTINUE
772 END IF
773*
774 IF( ilvr ) THEN
775 CALL cggbak( balanc, 'R', n, ilo, ihi, lscale, rscale, n,
776 $ vr,
777 $ ldvr, ierr )
778 DO 80 jc = 1, n
779 temp = zero
780 DO 60 jr = 1, n
781 temp = max( temp, abs1( vr( jr, jc ) ) )
782 60 CONTINUE
783 IF( temp.LT.smlnum )
784 $ GO TO 80
785 temp = one / temp
786 DO 70 jr = 1, n
787 vr( jr, jc ) = vr( jr, jc )*temp
788 70 CONTINUE
789 80 CONTINUE
790 END IF
791*
792* Undo scaling if necessary
793*
794 90 CONTINUE
795*
796 IF( ilascl )
797 $ CALL clascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
798*
799 IF( ilbscl )
800 $ CALL clascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
801*
802 work( 1 ) = sroundup_lwork(maxwrk)
803 RETURN
804*
805* End of CGGEVX
806*
807 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgeqrf(m, n, a, lda, tau, work, lwork, info)
CGEQRF
Definition cgeqrf.f:144
subroutine cggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
CGGBAK
Definition cggbak.f:147
subroutine cggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
CGGBAL
Definition cggbal.f:175
subroutine cggevx(balanc, jobvl, jobvr, sense, n, a, lda, b, ldb, alpha, beta, vl, ldvl, vr, ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, rconde, rcondv, work, lwork, rwork, iwork, bwork, info)
CGGEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
Definition cggevx.f:373
subroutine cgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
CGGHRD
Definition cgghrd.f:203
subroutine chgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
CHGEQZ
Definition chgeqz.f:283
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
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
subroutine ctgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
CTGEVC
Definition ctgevc.f:217
subroutine ctgsna(job, howmny, select, n, a, lda, b, ldb, vl, ldvl, vr, ldvr, s, dif, mm, m, work, lwork, iwork, info)
CTGSNA
Definition ctgsna.f:309
subroutine cungqr(m, n, k, a, lda, tau, work, lwork, info)
CUNGQR
Definition cungqr.f:126
subroutine cunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
CUNMQR
Definition cunmqr.f:166