LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zgeev()

subroutine zgeev ( character jobvl,
character jobvr,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( * ) w,
complex*16, dimension( ldvl, * ) vl,
integer ldvl,
complex*16, dimension( ldvr, * ) vr,
integer ldvr,
complex*16, dimension( * ) work,
integer lwork,
double precision, dimension( * ) rwork,
integer info )

ZGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices

Download ZGEEV + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> ZGEEV computes for an N-by-N complex nonsymmetric matrix A, the
!> eigenvalues and, optionally, the left and/or right eigenvectors.
!>
!> The right eigenvector v(j) of A satisfies
!>                  A * v(j) = lambda(j) * v(j)
!> where lambda(j) is its eigenvalue.
!> The left eigenvector u(j) of A satisfies
!>               u(j)**H * A = lambda(j) * u(j)**H
!> where u(j)**H denotes the conjugate transpose of u(j).
!>
!> The computed eigenvectors are normalized to have Euclidean norm
!> equal to 1 and largest component real.
!> 
Parameters
[in]JOBVL
!>          JOBVL is CHARACTER*1
!>          = 'N': left eigenvectors of A are not computed;
!>          = 'V': left eigenvectors of are computed.
!> 
[in]JOBVR
!>          JOBVR is CHARACTER*1
!>          = 'N': right eigenvectors of A are not computed;
!>          = 'V': right eigenvectors of A are computed.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A. N >= 0.
!> 
[in,out]A
!>          A is COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the N-by-N matrix A.
!>          On exit, A has been overwritten.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]W
!>          W is COMPLEX*16 array, dimension (N)
!>          W contains the computed eigenvalues.
!> 
[out]VL
!>          VL is COMPLEX*16 array, dimension (LDVL,N)
!>          If JOBVL = 'V', the left eigenvectors u(j) are stored one
!>          after another in the columns of VL, in the same order
!>          as their eigenvalues.
!>          If JOBVL = 'N', VL is not referenced.
!>          u(j) = VL(:,j), the j-th column of VL.
!> 
[in]LDVL
!>          LDVL is INTEGER
!>          The leading dimension of the array VL.  LDVL >= 1; if
!>          JOBVL = 'V', LDVL >= N.
!> 
[out]VR
!>          VR is COMPLEX*16 array, dimension (LDVR,N)
!>          If JOBVR = 'V', the right eigenvectors v(j) are stored one
!>          after another in the columns of VR, in the same order
!>          as their eigenvalues.
!>          If JOBVR = 'N', VR is not referenced.
!>          v(j) = VR(:,j), the j-th column of VR.
!> 
[in]LDVR
!>          LDVR is INTEGER
!>          The leading dimension of the array VR.  LDVR >= 1; if
!>          JOBVR = 'V', LDVR >= N.
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= max(1,2*N).
!>          For good performance, LWORK must generally be larger.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (2*N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!>          > 0:  if INFO = i, the QR algorithm failed to compute all the
!>                eigenvalues, and no eigenvectors have been computed;
!>                elements i+1:N of W contain eigenvalues which have
!>                converged.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 176 of file zgeev.f.

179 implicit none
180*
181* -- LAPACK driver routine --
182* -- LAPACK is a software package provided by Univ. of Tennessee, --
183* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
184*
185* .. Scalar Arguments ..
186 CHARACTER JOBVL, JOBVR
187 INTEGER INFO, LDA, LDVL, LDVR, LWORK, N
188* ..
189* .. Array Arguments ..
190 DOUBLE PRECISION RWORK( * )
191 COMPLEX*16 A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ),
192 $ W( * ), WORK( * )
193* ..
194*
195* =====================================================================
196*
197* .. Parameters ..
198 DOUBLE PRECISION ZERO, ONE
199 parameter( zero = 0.0d0, one = 1.0d0 )
200* ..
201* .. Local Scalars ..
202 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR
203 CHARACTER SIDE
204 INTEGER HSWORK, I, IBAL, IERR, IHI, ILO, IRWORK, ITAU,
205 $ IWRK, K, LWORK_TREVC, MAXWRK, MINWRK, NOUT
206 DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, SCL, SMLNUM
207 COMPLEX*16 TMP
208* ..
209* .. Local Arrays ..
210 LOGICAL SELECT( 1 )
211 DOUBLE PRECISION DUM( 1 )
212* ..
213* .. External Subroutines ..
214 EXTERNAL xerbla, zdscal, zgebak, zgebal, zgehrd,
215 $ zhseqr,
217* ..
218* .. External Functions ..
219 LOGICAL LSAME
220 INTEGER IDAMAX, ILAENV
221 DOUBLE PRECISION DLAMCH, DZNRM2, ZLANGE
222 EXTERNAL lsame, idamax, ilaenv, dlamch, dznrm2,
223 $ zlange
224* ..
225* .. Intrinsic Functions ..
226 INTRINSIC dble, dcmplx, conjg, aimag, max, sqrt
227* ..
228* .. Executable Statements ..
229*
230* Test the input arguments
231*
232 info = 0
233 lquery = ( lwork.EQ.-1 )
234 wantvl = lsame( jobvl, 'V' )
235 wantvr = lsame( jobvr, 'V' )
236 IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
237 info = -1
238 ELSE IF( ( .NOT.wantvr ) .AND.
239 $ ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
240 info = -2
241 ELSE IF( n.LT.0 ) THEN
242 info = -3
243 ELSE IF( lda.LT.max( 1, n ) ) THEN
244 info = -5
245 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
246 info = -8
247 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
248 info = -10
249 END IF
250*
251* Compute workspace
252* (Note: Comments in the code beginning "Workspace:" describe the
253* minimal amount of workspace needed at that point in the code,
254* as well as the preferred amount for good performance.
255* CWorkspace refers to complex workspace, and RWorkspace to real
256* workspace. NB refers to the optimal block size for the
257* immediately following subroutine, as returned by ILAENV.
258* HSWORK refers to the workspace preferred by ZHSEQR, as
259* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
260* the worst case.)
261*
262 IF( info.EQ.0 ) THEN
263 IF( n.EQ.0 ) THEN
264 minwrk = 1
265 maxwrk = 1
266 ELSE
267 maxwrk = n + n*ilaenv( 1, 'ZGEHRD', ' ', n, 1, n, 0 )
268 minwrk = 2*n
269 IF( wantvl ) THEN
270 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
271 $ 'ZUNGHR',
272 $ ' ', n, 1, n, -1 ) )
273 CALL ztrevc3( 'L', 'B', SELECT, n, a, lda,
274 $ vl, ldvl, vr, ldvr,
275 $ n, nout, work, -1, rwork, -1, ierr )
276 lwork_trevc = int( work(1) )
277 maxwrk = max( maxwrk, n + lwork_trevc )
278 CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vl, ldvl,
279 $ work, -1, info )
280 ELSE IF( wantvr ) THEN
281 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
282 $ 'ZUNGHR',
283 $ ' ', n, 1, n, -1 ) )
284 CALL ztrevc3( 'R', 'B', SELECT, n, a, lda,
285 $ vl, ldvl, vr, ldvr,
286 $ n, nout, work, -1, rwork, -1, ierr )
287 lwork_trevc = int( work(1) )
288 maxwrk = max( maxwrk, n + lwork_trevc )
289 CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vr, ldvr,
290 $ work, -1, info )
291 ELSE
292 CALL zhseqr( 'E', 'N', n, 1, n, a, lda, w, vr, ldvr,
293 $ work, -1, info )
294 END IF
295 hswork = int( work(1) )
296 maxwrk = max( maxwrk, hswork, minwrk )
297 END IF
298 work( 1 ) = maxwrk
299*
300 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
301 info = -12
302 END IF
303 END IF
304*
305 IF( info.NE.0 ) THEN
306 CALL xerbla( 'ZGEEV ', -info )
307 RETURN
308 ELSE IF( lquery ) THEN
309 RETURN
310 END IF
311*
312* Quick return if possible
313*
314 IF( n.EQ.0 )
315 $ RETURN
316*
317* Get machine constants
318*
319 eps = dlamch( 'P' )
320 smlnum = dlamch( 'S' )
321 bignum = one / smlnum
322 smlnum = sqrt( smlnum ) / eps
323 bignum = one / smlnum
324*
325* Scale A if max element outside range [SMLNUM,BIGNUM]
326*
327 anrm = zlange( 'M', n, n, a, lda, dum )
328 scalea = .false.
329 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
330 scalea = .true.
331 cscale = smlnum
332 ELSE IF( anrm.GT.bignum ) THEN
333 scalea = .true.
334 cscale = bignum
335 END IF
336 IF( scalea )
337 $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
338*
339* Balance the matrix
340* (CWorkspace: none)
341* (RWorkspace: need N)
342*
343 ibal = 1
344 CALL zgebal( 'B', n, a, lda, ilo, ihi, rwork( ibal ), ierr )
345*
346* Reduce to upper Hessenberg form
347* (CWorkspace: need 2*N, prefer N+N*NB)
348* (RWorkspace: none)
349*
350 itau = 1
351 iwrk = itau + n
352 CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
353 $ lwork-iwrk+1, ierr )
354*
355 IF( wantvl ) THEN
356*
357* Want left eigenvectors
358* Copy Householder vectors to VL
359*
360 side = 'L'
361 CALL zlacpy( 'L', n, n, a, lda, vl, ldvl )
362*
363* Generate unitary matrix in VL
364* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
365* (RWorkspace: none)
366*
367 CALL zunghr( n, ilo, ihi, vl, ldvl, work( itau ),
368 $ work( iwrk ),
369 $ lwork-iwrk+1, ierr )
370*
371* Perform QR iteration, accumulating Schur vectors in VL
372* (CWorkspace: need 1, prefer HSWORK (see comments) )
373* (RWorkspace: none)
374*
375 iwrk = itau
376 CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vl, ldvl,
377 $ work( iwrk ), lwork-iwrk+1, info )
378*
379 IF( wantvr ) THEN
380*
381* Want left and right eigenvectors
382* Copy Schur vectors to VR
383*
384 side = 'B'
385 CALL zlacpy( 'F', n, n, vl, ldvl, vr, ldvr )
386 END IF
387*
388 ELSE IF( wantvr ) THEN
389*
390* Want right eigenvectors
391* Copy Householder vectors to VR
392*
393 side = 'R'
394 CALL zlacpy( 'L', n, n, a, lda, vr, ldvr )
395*
396* Generate unitary matrix in VR
397* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
398* (RWorkspace: none)
399*
400 CALL zunghr( n, ilo, ihi, vr, ldvr, work( itau ),
401 $ work( iwrk ),
402 $ lwork-iwrk+1, ierr )
403*
404* Perform QR iteration, accumulating Schur vectors in VR
405* (CWorkspace: need 1, prefer HSWORK (see comments) )
406* (RWorkspace: none)
407*
408 iwrk = itau
409 CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vr, ldvr,
410 $ work( iwrk ), lwork-iwrk+1, info )
411*
412 ELSE
413*
414* Compute eigenvalues only
415* (CWorkspace: need 1, prefer HSWORK (see comments) )
416* (RWorkspace: none)
417*
418 iwrk = itau
419 CALL zhseqr( 'E', 'N', n, ilo, ihi, a, lda, w, vr, ldvr,
420 $ work( iwrk ), lwork-iwrk+1, info )
421 END IF
422*
423* If INFO .NE. 0 from ZHSEQR, then quit
424*
425 IF( info.NE.0 )
426 $ GO TO 50
427*
428 IF( wantvl .OR. wantvr ) THEN
429*
430* Compute left and/or right eigenvectors
431* (CWorkspace: need 2*N, prefer N + 2*N*NB)
432* (RWorkspace: need 2*N)
433*
434 irwork = ibal + n
435 CALL ztrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr,
436 $ ldvr,
437 $ n, nout, work( iwrk ), lwork-iwrk+1,
438 $ rwork( irwork ), n, ierr )
439 END IF
440*
441 IF( wantvl ) THEN
442*
443* Undo balancing of left eigenvectors
444* (CWorkspace: none)
445* (RWorkspace: need N)
446*
447 CALL zgebak( 'B', 'L', n, ilo, ihi, rwork( ibal ), n, vl,
448 $ ldvl,
449 $ ierr )
450*
451* Normalize left eigenvectors and make largest component real
452*
453 DO 20 i = 1, n
454 scl = one / dznrm2( n, vl( 1, i ), 1 )
455 CALL zdscal( n, scl, vl( 1, i ), 1 )
456 DO 10 k = 1, n
457 rwork( irwork+k-1 ) = dble( vl( k, i ) )**2 +
458 $ aimag( vl( k, i ) )**2
459 10 CONTINUE
460 k = idamax( n, rwork( irwork ), 1 )
461 tmp = conjg( vl( k, i ) ) / sqrt( rwork( irwork+k-1 ) )
462 CALL zscal( n, tmp, vl( 1, i ), 1 )
463 vl( k, i ) = dcmplx( dble( vl( k, i ) ), zero )
464 20 CONTINUE
465 END IF
466*
467 IF( wantvr ) THEN
468*
469* Undo balancing of right eigenvectors
470* (CWorkspace: none)
471* (RWorkspace: need N)
472*
473 CALL zgebak( 'B', 'R', n, ilo, ihi, rwork( ibal ), n, vr,
474 $ ldvr,
475 $ ierr )
476*
477* Normalize right eigenvectors and make largest component real
478*
479 DO 40 i = 1, n
480 scl = one / dznrm2( n, vr( 1, i ), 1 )
481 CALL zdscal( n, scl, vr( 1, i ), 1 )
482 DO 30 k = 1, n
483 rwork( irwork+k-1 ) = dble( vr( k, i ) )**2 +
484 $ aimag( vr( k, i ) )**2
485 30 CONTINUE
486 k = idamax( n, rwork( irwork ), 1 )
487 tmp = conjg( vr( k, i ) ) / sqrt( rwork( irwork+k-1 ) )
488 CALL zscal( n, tmp, vr( 1, i ), 1 )
489 vr( k, i ) = dcmplx( dble( vr( k, i ) ), zero )
490 40 CONTINUE
491 END IF
492*
493* Undo scaling if necessary
494*
495 50 CONTINUE
496 IF( scalea ) THEN
497 CALL zlascl( 'G', 0, 0, cscale, anrm, n-info, 1,
498 $ w( info+1 ),
499 $ max( n-info, 1 ), ierr )
500 IF( info.GT.0 ) THEN
501 CALL zlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, w, n,
502 $ ierr )
503 END IF
504 END IF
505*
506 work( 1 ) = maxwrk
507 RETURN
508*
509* End of ZGEEV
510*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
ZGEBAK
Definition zgebak.f:129
subroutine zgebal(job, n, a, lda, ilo, ihi, scale, info)
ZGEBAL
Definition zgebal.f:163
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
Definition zgehrd.f:166
subroutine zhseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
ZHSEQR
Definition zhseqr.f:297
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:113
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function dznrm2(n, x, incx)
DZNRM2
Definition dznrm2.f90:90
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:243
subroutine zunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZUNGHR
Definition zunghr.f:125
Here is the call graph for this function:
Here is the caller graph for this function: