LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zggev.f
Go to the documentation of this file.
1*> \brief <b> ZGGEV 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 ZGGEV + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zggev.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zggev.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zggev.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGGEV( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHA, BETA,
20* VL, LDVL, VR, LDVR, WORK, LWORK, RWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER JOBVL, JOBVR
24* INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N
25* ..
26* .. Array Arguments ..
27* DOUBLE PRECISION RWORK( * )
28* COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
29* $ BETA( * ), VL( LDVL, * ), VR( LDVR, * ),
30* $ WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> ZGGEV computes for a pair of N-by-N complex nonsymmetric matrices
40*> (A,B), the generalized eigenvalues, and optionally, the left and/or
41*> right generalized eigenvectors.
42*>
43*> A generalized eigenvalue for a pair of matrices (A,B) is a scalar
44*> lambda or a ratio alpha/beta = lambda, such that A - lambda*B is
45*> singular. It is usually represented as the pair (alpha,beta), as
46*> there is a reasonable interpretation for beta=0, and even for both
47*> being zero.
48*>
49*> The right generalized eigenvector v(j) corresponding to the
50*> generalized eigenvalue lambda(j) of (A,B) satisfies
51*>
52*> A * v(j) = lambda(j) * B * v(j).
53*>
54*> The left generalized eigenvector u(j) corresponding to the
55*> generalized eigenvalues lambda(j) of (A,B) satisfies
56*>
57*> u(j)**H * A = lambda(j) * u(j)**H * B
58*>
59*> where u(j)**H is the conjugate-transpose of u(j).
60*> \endverbatim
61*
62* Arguments:
63* ==========
64*
65*> \param[in] JOBVL
66*> \verbatim
67*> JOBVL is CHARACTER*1
68*> = 'N': do not compute the left generalized eigenvectors;
69*> = 'V': compute the left generalized eigenvectors.
70*> \endverbatim
71*>
72*> \param[in] JOBVR
73*> \verbatim
74*> JOBVR is CHARACTER*1
75*> = 'N': do not compute the right generalized eigenvectors;
76*> = 'V': compute the right generalized eigenvectors.
77*> \endverbatim
78*>
79*> \param[in] N
80*> \verbatim
81*> N is INTEGER
82*> The order of the matrices A, B, VL, and VR. N >= 0.
83*> \endverbatim
84*>
85*> \param[in,out] A
86*> \verbatim
87*> A is COMPLEX*16 array, dimension (LDA, N)
88*> On entry, the matrix A in the pair (A,B).
89*> On exit, A has been overwritten.
90*> \endverbatim
91*>
92*> \param[in] LDA
93*> \verbatim
94*> LDA is INTEGER
95*> The leading dimension of A. LDA >= max(1,N).
96*> \endverbatim
97*>
98*> \param[in,out] B
99*> \verbatim
100*> B is COMPLEX*16 array, dimension (LDB, N)
101*> On entry, the matrix B in the pair (A,B).
102*> On exit, B has been overwritten.
103*> \endverbatim
104*>
105*> \param[in] LDB
106*> \verbatim
107*> LDB is INTEGER
108*> The leading dimension of B. LDB >= max(1,N).
109*> \endverbatim
110*>
111*> \param[out] ALPHA
112*> \verbatim
113*> ALPHA is COMPLEX*16 array, dimension (N)
114*> \endverbatim
115*>
116*> \param[out] BETA
117*> \verbatim
118*> BETA is COMPLEX*16 array, dimension (N)
119*> On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the
120*> generalized eigenvalues.
121*>
122*> Note: the quotients ALPHA(j)/BETA(j) may easily over- or
123*> underflow, and BETA(j) may even be zero. Thus, the user
124*> should avoid naively computing the ratio alpha/beta.
125*> However, ALPHA will be always less than and usually
126*> comparable with norm(A) in magnitude, and BETA always less
127*> than and usually comparable with norm(B).
128*> \endverbatim
129*>
130*> \param[out] VL
131*> \verbatim
132*> VL is COMPLEX*16 array, dimension (LDVL,N)
133*> If JOBVL = 'V', the left generalized eigenvectors u(j) are
134*> stored one after another in the columns of VL, in the same
135*> order as their eigenvalues.
136*> Each eigenvector is scaled so the largest component has
137*> abs(real part) + abs(imag. part) = 1.
138*> Not referenced if JOBVL = 'N'.
139*> \endverbatim
140*>
141*> \param[in] LDVL
142*> \verbatim
143*> LDVL is INTEGER
144*> The leading dimension of the matrix VL. LDVL >= 1, and
145*> if JOBVL = 'V', LDVL >= N.
146*> \endverbatim
147*>
148*> \param[out] VR
149*> \verbatim
150*> VR is COMPLEX*16 array, dimension (LDVR,N)
151*> If JOBVR = 'V', the right generalized eigenvectors v(j) are
152*> stored one after another in the columns of VR, in the same
153*> order as their eigenvalues.
154*> Each eigenvector is scaled so the largest component has
155*> abs(real part) + abs(imag. part) = 1.
156*> Not referenced if JOBVR = 'N'.
157*> \endverbatim
158*>
159*> \param[in] LDVR
160*> \verbatim
161*> LDVR is INTEGER
162*> The leading dimension of the matrix VR. LDVR >= 1, and
163*> if JOBVR = 'V', LDVR >= N.
164*> \endverbatim
165*>
166*> \param[out] WORK
167*> \verbatim
168*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
169*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
170*> \endverbatim
171*>
172*> \param[in] LWORK
173*> \verbatim
174*> LWORK is INTEGER
175*> The dimension of the array WORK. LWORK >= max(1,2*N).
176*> For good performance, LWORK must generally be larger.
177*>
178*> If LWORK = -1, then a workspace query is assumed; the routine
179*> only calculates the optimal size of the WORK array, returns
180*> this value as the first entry of the WORK array, and no error
181*> message related to LWORK is issued by XERBLA.
182*> \endverbatim
183*>
184*> \param[out] RWORK
185*> \verbatim
186*> RWORK is DOUBLE PRECISION array, dimension (8*N)
187*> \endverbatim
188*>
189*> \param[out] INFO
190*> \verbatim
191*> INFO is INTEGER
192*> = 0: successful exit
193*> < 0: if INFO = -i, the i-th argument had an illegal value.
194*> =1,...,N:
195*> The QZ iteration failed. No eigenvectors have been
196*> calculated, but ALPHA(j) and BETA(j) should be
197*> correct for j=INFO+1,...,N.
198*> > N: =N+1: other then QZ iteration failed in ZHGEQZ,
199*> =N+2: error return from ZTGEVC.
200*> \endverbatim
201*
202* Authors:
203* ========
204*
205*> \author Univ. of Tennessee
206*> \author Univ. of California Berkeley
207*> \author Univ. of Colorado Denver
208*> \author NAG Ltd.
209*
210*> \ingroup ggev
211*
212* =====================================================================
213 SUBROUTINE zggev( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHA, BETA,
214 $ VL, LDVL, VR, LDVR, WORK, LWORK, RWORK, INFO )
215*
216* -- LAPACK driver routine --
217* -- LAPACK is a software package provided by Univ. of Tennessee, --
218* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
219*
220* .. Scalar Arguments ..
221 CHARACTER JOBVL, JOBVR
222 INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N
223* ..
224* .. Array Arguments ..
225 DOUBLE PRECISION RWORK( * )
226 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
227 $ beta( * ), vl( ldvl, * ), vr( ldvr, * ),
228 $ work( * )
229* ..
230*
231* =====================================================================
232*
233* .. Parameters ..
234 DOUBLE PRECISION ZERO, ONE
235 parameter( zero = 0.0d0, one = 1.0d0 )
236 COMPLEX*16 CZERO, CONE
237 parameter( czero = ( 0.0d0, 0.0d0 ),
238 $ cone = ( 1.0d0, 0.0d0 ) )
239* ..
240* .. Local Scalars ..
241 LOGICAL ILASCL, ILBSCL, ILV, ILVL, ILVR, LQUERY
242 CHARACTER CHTEMP
243 INTEGER ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT, ILO,
244 $ in, iright, irows, irwrk, itau, iwrk, jc, jr,
245 $ lwkmin, lwkopt
246 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
247 $ smlnum, temp
248 COMPLEX*16 X
249* ..
250* .. Local Arrays ..
251 LOGICAL LDUMMA( 1 )
252* ..
253* .. External Subroutines ..
254 EXTERNAL xerbla, zgeqrf, zggbak, zggbal, zgghrd,
255 $ zhgeqz,
257* ..
258* .. External Functions ..
259 LOGICAL LSAME
260 INTEGER ILAENV
261 DOUBLE PRECISION DLAMCH, ZLANGE
262 EXTERNAL lsame, ilaenv, dlamch, zlange
263* ..
264* .. Intrinsic Functions ..
265 INTRINSIC abs, dble, dimag, max, sqrt
266* ..
267* .. Statement Functions ..
268 DOUBLE PRECISION ABS1
269* ..
270* .. Statement Function definitions ..
271 abs1( x ) = abs( dble( x ) ) + abs( dimag( x ) )
272* ..
273* .. Executable Statements ..
274*
275* Decode the input arguments
276*
277 IF( lsame( jobvl, 'N' ) ) THEN
278 ijobvl = 1
279 ilvl = .false.
280 ELSE IF( lsame( jobvl, 'V' ) ) THEN
281 ijobvl = 2
282 ilvl = .true.
283 ELSE
284 ijobvl = -1
285 ilvl = .false.
286 END IF
287*
288 IF( lsame( jobvr, 'N' ) ) THEN
289 ijobvr = 1
290 ilvr = .false.
291 ELSE IF( lsame( jobvr, 'V' ) ) THEN
292 ijobvr = 2
293 ilvr = .true.
294 ELSE
295 ijobvr = -1
296 ilvr = .false.
297 END IF
298 ilv = ilvl .OR. ilvr
299*
300* Test the input arguments
301*
302 info = 0
303 lquery = ( lwork.EQ.-1 )
304 IF( ijobvl.LE.0 ) THEN
305 info = -1
306 ELSE IF( ijobvr.LE.0 ) THEN
307 info = -2
308 ELSE IF( n.LT.0 ) THEN
309 info = -3
310 ELSE IF( lda.LT.max( 1, n ) ) THEN
311 info = -5
312 ELSE IF( ldb.LT.max( 1, n ) ) THEN
313 info = -7
314 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
315 info = -11
316 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
317 info = -13
318 END IF
319*
320* Compute workspace
321* (Note: Comments in the code beginning "Workspace:" describe the
322* minimal amount of workspace needed at that point in the code,
323* as well as the preferred amount for good performance.
324* NB refers to the optimal block size for the immediately
325* following subroutine, as returned by ILAENV. The workspace is
326* computed assuming ILO = 1 and IHI = N, the worst case.)
327*
328 IF( info.EQ.0 ) THEN
329 lwkmin = max( 1, 2*n )
330 lwkopt = max( 1, n + n*ilaenv( 1, 'ZGEQRF', ' ', n, 1, n,
331 $ 0 ) )
332 lwkopt = max( lwkopt, n +
333 $ n*ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, 0 ) )
334 IF( ilvl ) THEN
335 lwkopt = max( lwkopt, n +
336 $ n*ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, -1 ) )
337 END IF
338 work( 1 ) = lwkopt
339*
340 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
341 $ info = -15
342 END IF
343*
344 IF( info.NE.0 ) THEN
345 CALL xerbla( 'ZGGEV ', -info )
346 RETURN
347 ELSE IF( lquery ) THEN
348 RETURN
349 END IF
350*
351* Quick return if possible
352*
353 IF( n.EQ.0 )
354 $ RETURN
355*
356* Get machine constants
357*
358 eps = dlamch( 'E' )*dlamch( 'B' )
359 smlnum = dlamch( 'S' )
360 bignum = one / smlnum
361 smlnum = sqrt( smlnum ) / eps
362 bignum = one / smlnum
363*
364* Scale A if max element outside range [SMLNUM,BIGNUM]
365*
366 anrm = zlange( 'M', n, n, a, lda, rwork )
367 ilascl = .false.
368 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
369 anrmto = smlnum
370 ilascl = .true.
371 ELSE IF( anrm.GT.bignum ) THEN
372 anrmto = bignum
373 ilascl = .true.
374 END IF
375 IF( ilascl )
376 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
377*
378* Scale B if max element outside range [SMLNUM,BIGNUM]
379*
380 bnrm = zlange( 'M', n, n, b, ldb, rwork )
381 ilbscl = .false.
382 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
383 bnrmto = smlnum
384 ilbscl = .true.
385 ELSE IF( bnrm.GT.bignum ) THEN
386 bnrmto = bignum
387 ilbscl = .true.
388 END IF
389 IF( ilbscl )
390 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
391*
392* Permute the matrices A, B to isolate eigenvalues if possible
393* (Real Workspace: need 6*N)
394*
395 ileft = 1
396 iright = n + 1
397 irwrk = iright + n
398 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
399 $ rwork( iright ), rwork( irwrk ), ierr )
400*
401* Reduce B to triangular form (QR decomposition of B)
402* (Complex Workspace: need N, prefer N*NB)
403*
404 irows = ihi + 1 - ilo
405 IF( ilv ) THEN
406 icols = n + 1 - ilo
407 ELSE
408 icols = irows
409 END IF
410 itau = 1
411 iwrk = itau + irows
412 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
413 $ work( iwrk ), lwork+1-iwrk, ierr )
414*
415* Apply the orthogonal transformation to matrix A
416* (Complex Workspace: need N, prefer N*NB)
417*
418 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
419 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
420 $ lwork+1-iwrk, ierr )
421*
422* Initialize VL
423* (Complex Workspace: need N, prefer N*NB)
424*
425 IF( ilvl ) THEN
426 CALL zlaset( 'Full', n, n, czero, cone, vl, ldvl )
427 IF( irows.GT.1 ) THEN
428 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
429 $ vl( ilo+1, ilo ), ldvl )
430 END IF
431 CALL zungqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
432 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
433 END IF
434*
435* Initialize VR
436*
437 IF( ilvr )
438 $ CALL zlaset( 'Full', n, n, czero, cone, vr, ldvr )
439*
440* Reduce to generalized Hessenberg form
441*
442 IF( ilv ) THEN
443*
444* Eigenvectors requested -- work on whole matrix.
445*
446 CALL zgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
447 $ ldvl, vr, ldvr, ierr )
448 ELSE
449 CALL zgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
450 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
451 END IF
452*
453* Perform QZ algorithm (Compute eigenvalues, and optionally, the
454* Schur form and Schur vectors)
455* (Complex Workspace: need N)
456* (Real Workspace: need N)
457*
458 iwrk = itau
459 IF( ilv ) THEN
460 chtemp = 'S'
461 ELSE
462 chtemp = 'E'
463 END IF
464 CALL zhgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
465 $ alpha, beta, vl, ldvl, vr, ldvr, work( iwrk ),
466 $ lwork+1-iwrk, rwork( irwrk ), ierr )
467 IF( ierr.NE.0 ) THEN
468 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
469 info = ierr
470 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
471 info = ierr - n
472 ELSE
473 info = n + 1
474 END IF
475 GO TO 70
476 END IF
477*
478* Compute Eigenvectors
479* (Real Workspace: need 2*N)
480* (Complex Workspace: need 2*N)
481*
482 IF( ilv ) THEN
483 IF( ilvl ) THEN
484 IF( ilvr ) THEN
485 chtemp = 'B'
486 ELSE
487 chtemp = 'L'
488 END IF
489 ELSE
490 chtemp = 'R'
491 END IF
492*
493 CALL ztgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl,
494 $ ldvl,
495 $ vr, ldvr, n, in, work( iwrk ), rwork( irwrk ),
496 $ ierr )
497 IF( ierr.NE.0 ) THEN
498 info = n + 2
499 GO TO 70
500 END IF
501*
502* Undo balancing on VL and VR and normalization
503* (Workspace: none needed)
504*
505 IF( ilvl ) THEN
506 CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
507 $ rwork( iright ), n, vl, ldvl, ierr )
508 DO 30 jc = 1, n
509 temp = zero
510 DO 10 jr = 1, n
511 temp = max( temp, abs1( vl( jr, jc ) ) )
512 10 CONTINUE
513 IF( temp.LT.smlnum )
514 $ GO TO 30
515 temp = one / temp
516 DO 20 jr = 1, n
517 vl( jr, jc ) = vl( jr, jc )*temp
518 20 CONTINUE
519 30 CONTINUE
520 END IF
521 IF( ilvr ) THEN
522 CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
523 $ rwork( iright ), n, vr, ldvr, ierr )
524 DO 60 jc = 1, n
525 temp = zero
526 DO 40 jr = 1, n
527 temp = max( temp, abs1( vr( jr, jc ) ) )
528 40 CONTINUE
529 IF( temp.LT.smlnum )
530 $ GO TO 60
531 temp = one / temp
532 DO 50 jr = 1, n
533 vr( jr, jc ) = vr( jr, jc )*temp
534 50 CONTINUE
535 60 CONTINUE
536 END IF
537 END IF
538*
539* Undo scaling if necessary
540*
541 70 CONTINUE
542*
543 IF( ilascl )
544 $ CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
545*
546 IF( ilbscl )
547 $ CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
548*
549 work( 1 ) = lwkopt
550 RETURN
551*
552* End of ZGGEV
553*
554 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqrf(m, n, a, lda, tau, work, lwork, info)
ZGEQRF
Definition zgeqrf.f:144
subroutine zggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
ZGGBAK
Definition zggbak.f:147
subroutine zggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
ZGGBAL
Definition zggbal.f:175
subroutine zggev(jobvl, jobvr, n, a, lda, b, ldb, alpha, beta, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
ZGGEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
Definition zggev.f:215
subroutine zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:203
subroutine zhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
ZHGEQZ
Definition zhgeqz.f:283
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
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:142
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:104
subroutine ztgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
ZTGEVC
Definition ztgevc.f:217
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:126
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:165