LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 and i+1:N of W contain eigenvalues which have
                converged.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
June 2016

Definition at line 181 of file zgeev.f.

181  implicit none
182 *
183 * -- LAPACK driver routine (version 3.6.1) --
184 * -- LAPACK is a software package provided by Univ. of Tennessee, --
185 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
186 * June 2016
187 *
188 * .. Scalar Arguments ..
189  CHARACTER jobvl, jobvr
190  INTEGER info, lda, ldvl, ldvr, lwork, n
191 * ..
192 * .. Array Arguments ..
193  DOUBLE PRECISION rwork( * )
194  COMPLEX*16 a( lda, * ), vl( ldvl, * ), vr( ldvr, * ),
195  $ w( * ), work( * )
196 * ..
197 *
198 * =====================================================================
199 *
200 * .. Parameters ..
201  DOUBLE PRECISION zero, one
202  parameter ( zero = 0.0d0, one = 1.0d0 )
203 * ..
204 * .. Local Scalars ..
205  LOGICAL lquery, scalea, wantvl, wantvr
206  CHARACTER side
207  INTEGER hswork, i, ibal, ierr, ihi, ilo, irwork, itau,
208  $ iwrk, k, lwork_trevc, maxwrk, minwrk, nout
209  DOUBLE PRECISION anrm, bignum, cscale, eps, scl, smlnum
210  COMPLEX*16 tmp
211 * ..
212 * .. Local Arrays ..
213  LOGICAL select( 1 )
214  DOUBLE PRECISION dum( 1 )
215 * ..
216 * .. External Subroutines ..
217  EXTERNAL dlabad, xerbla, zdscal, zgebak, zgebal, zgehrd,
219 * ..
220 * .. External Functions ..
221  LOGICAL lsame
222  INTEGER idamax, ilaenv
223  DOUBLE PRECISION dlamch, dznrm2, zlange
224  EXTERNAL lsame, idamax, ilaenv, dlamch, dznrm2, zlange
225 * ..
226 * .. Intrinsic Functions ..
227  INTRINSIC dble, dcmplx, conjg, aimag, max, sqrt
228 * ..
229 * .. Executable Statements ..
230 *
231 * Test the input arguments
232 *
233  info = 0
234  lquery = ( lwork.EQ.-1 )
235  wantvl = lsame( jobvl, 'V' )
236  wantvr = lsame( jobvr, 'V' )
237  IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
238  info = -1
239  ELSE IF( ( .NOT.wantvr ) .AND. ( .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, 'ZUNGHR',
271  $ ' ', n, 1, n, -1 ) )
272  CALL ztrevc3( 'L', 'B', SELECT, n, a, lda,
273  $ vl, ldvl, vr, ldvr,
274  $ n, nout, work, -1, rwork, -1, ierr )
275  lwork_trevc = int( work(1) )
276  maxwrk = max( maxwrk, n + lwork_trevc )
277  CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vl, ldvl,
278  $ work, -1, info )
279  ELSE IF( wantvr ) THEN
280  maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1, 'ZUNGHR',
281  $ ' ', n, 1, n, -1 ) )
282  CALL ztrevc3( 'R', 'B', SELECT, n, a, lda,
283  $ vl, ldvl, vr, ldvr,
284  $ n, nout, work, -1, rwork, -1, ierr )
285  lwork_trevc = int( work(1) )
286  maxwrk = max( maxwrk, n + lwork_trevc )
287  CALL zhseqr( 'S', 'V', n, 1, n, a, lda, w, vr, ldvr,
288  $ work, -1, info )
289  ELSE
290  CALL zhseqr( 'E', 'N', n, 1, n, a, lda, w, vr, ldvr,
291  $ work, -1, info )
292  END IF
293  hswork = int( work(1) )
294  maxwrk = max( maxwrk, hswork, minwrk )
295  END IF
296  work( 1 ) = maxwrk
297 *
298  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
299  info = -12
300  END IF
301  END IF
302 *
303  IF( info.NE.0 ) THEN
304  CALL xerbla( 'ZGEEV ', -info )
305  RETURN
306  ELSE IF( lquery ) THEN
307  RETURN
308  END IF
309 *
310 * Quick return if possible
311 *
312  IF( n.EQ.0 )
313  $ RETURN
314 *
315 * Get machine constants
316 *
317  eps = dlamch( 'P' )
318  smlnum = dlamch( 'S' )
319  bignum = one / smlnum
320  CALL dlabad( smlnum, bignum )
321  smlnum = sqrt( smlnum ) / eps
322  bignum = one / smlnum
323 *
324 * Scale A if max element outside range [SMLNUM,BIGNUM]
325 *
326  anrm = zlange( 'M', n, n, a, lda, dum )
327  scalea = .false.
328  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
329  scalea = .true.
330  cscale = smlnum
331  ELSE IF( anrm.GT.bignum ) THEN
332  scalea = .true.
333  cscale = bignum
334  END IF
335  IF( scalea )
336  $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
337 *
338 * Balance the matrix
339 * (CWorkspace: none)
340 * (RWorkspace: need N)
341 *
342  ibal = 1
343  CALL zgebal( 'B', n, a, lda, ilo, ihi, rwork( ibal ), ierr )
344 *
345 * Reduce to upper Hessenberg form
346 * (CWorkspace: need 2*N, prefer N+N*NB)
347 * (RWorkspace: none)
348 *
349  itau = 1
350  iwrk = itau + n
351  CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
352  $ lwork-iwrk+1, ierr )
353 *
354  IF( wantvl ) THEN
355 *
356 * Want left eigenvectors
357 * Copy Householder vectors to VL
358 *
359  side = 'L'
360  CALL zlacpy( 'L', n, n, a, lda, vl, ldvl )
361 *
362 * Generate unitary matrix in VL
363 * (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
364 * (RWorkspace: none)
365 *
366  CALL zunghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
367  $ lwork-iwrk+1, ierr )
368 *
369 * Perform QR iteration, accumulating Schur vectors in VL
370 * (CWorkspace: need 1, prefer HSWORK (see comments) )
371 * (RWorkspace: none)
372 *
373  iwrk = itau
374  CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vl, ldvl,
375  $ work( iwrk ), lwork-iwrk+1, info )
376 *
377  IF( wantvr ) THEN
378 *
379 * Want left and right eigenvectors
380 * Copy Schur vectors to VR
381 *
382  side = 'B'
383  CALL zlacpy( 'F', n, n, vl, ldvl, vr, ldvr )
384  END IF
385 *
386  ELSE IF( wantvr ) THEN
387 *
388 * Want right eigenvectors
389 * Copy Householder vectors to VR
390 *
391  side = 'R'
392  CALL zlacpy( 'L', n, n, a, lda, vr, ldvr )
393 *
394 * Generate unitary matrix in VR
395 * (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
396 * (RWorkspace: none)
397 *
398  CALL zunghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
399  $ lwork-iwrk+1, ierr )
400 *
401 * Perform QR iteration, accumulating Schur vectors in VR
402 * (CWorkspace: need 1, prefer HSWORK (see comments) )
403 * (RWorkspace: none)
404 *
405  iwrk = itau
406  CALL zhseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vr, ldvr,
407  $ work( iwrk ), lwork-iwrk+1, info )
408 *
409  ELSE
410 *
411 * Compute eigenvalues only
412 * (CWorkspace: need 1, prefer HSWORK (see comments) )
413 * (RWorkspace: none)
414 *
415  iwrk = itau
416  CALL zhseqr( 'E', 'N', n, ilo, ihi, a, lda, w, vr, ldvr,
417  $ work( iwrk ), lwork-iwrk+1, info )
418  END IF
419 *
420 * If INFO .NE. 0 from ZHSEQR, then quit
421 *
422  IF( info.NE.0 )
423  $ GO TO 50
424 *
425  IF( wantvl .OR. wantvr ) THEN
426 *
427 * Compute left and/or right eigenvectors
428 * (CWorkspace: need 2*N, prefer N + 2*N*NB)
429 * (RWorkspace: need 2*N)
430 *
431  irwork = ibal + n
432  CALL ztrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
433  $ n, nout, work( iwrk ), lwork-iwrk+1,
434  $ rwork( irwork ), n, ierr )
435  END IF
436 *
437  IF( wantvl ) THEN
438 *
439 * Undo balancing of left eigenvectors
440 * (CWorkspace: none)
441 * (RWorkspace: need N)
442 *
443  CALL zgebak( 'B', 'L', n, ilo, ihi, rwork( ibal ), n, vl, ldvl,
444  $ ierr )
445 *
446 * Normalize left eigenvectors and make largest component real
447 *
448  DO 20 i = 1, n
449  scl = one / dznrm2( n, vl( 1, i ), 1 )
450  CALL zdscal( n, scl, vl( 1, i ), 1 )
451  DO 10 k = 1, n
452  rwork( irwork+k-1 ) = dble( vl( k, i ) )**2 +
453  $ aimag( vl( k, i ) )**2
454  10 CONTINUE
455  k = idamax( n, rwork( irwork ), 1 )
456  tmp = conjg( vl( k, i ) ) / sqrt( rwork( irwork+k-1 ) )
457  CALL zscal( n, tmp, vl( 1, i ), 1 )
458  vl( k, i ) = dcmplx( dble( vl( k, i ) ), zero )
459  20 CONTINUE
460  END IF
461 *
462  IF( wantvr ) THEN
463 *
464 * Undo balancing of right eigenvectors
465 * (CWorkspace: none)
466 * (RWorkspace: need N)
467 *
468  CALL zgebak( 'B', 'R', n, ilo, ihi, rwork( ibal ), n, vr, ldvr,
469  $ ierr )
470 *
471 * Normalize right eigenvectors and make largest component real
472 *
473  DO 40 i = 1, n
474  scl = one / dznrm2( n, vr( 1, i ), 1 )
475  CALL zdscal( n, scl, vr( 1, i ), 1 )
476  DO 30 k = 1, n
477  rwork( irwork+k-1 ) = dble( vr( k, i ) )**2 +
478  $ aimag( vr( k, i ) )**2
479  30 CONTINUE
480  k = idamax( n, rwork( irwork ), 1 )
481  tmp = conjg( vr( k, i ) ) / sqrt( rwork( irwork+k-1 ) )
482  CALL zscal( n, tmp, vr( 1, i ), 1 )
483  vr( k, i ) = dcmplx( dble( vr( k, i ) ), zero )
484  40 CONTINUE
485  END IF
486 *
487 * Undo scaling if necessary
488 *
489  50 CONTINUE
490  IF( scalea ) THEN
491  CALL zlascl( 'G', 0, 0, cscale, anrm, n-info, 1, w( info+1 ),
492  $ max( n-info, 1 ), ierr )
493  IF( info.GT.0 ) THEN
494  CALL zlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, w, n, ierr )
495  END IF
496  END IF
497 *
498  work( 1 ) = maxwrk
499  RETURN
500 *
501 * End of ZGEEV
502 *
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgebak(JOB, SIDE, N, ILO, IHI, SCALE, M, V, LDV, INFO)
ZGEBAK
Definition: zgebak.f:133
subroutine zgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZGEHRD
Definition: zgehrd.f:169
subroutine zunghr(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGHR
Definition: zunghr.f:128
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
double precision function dznrm2(N, X, INCX)
DZNRM2
Definition: dznrm2.f:56
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:117
subroutine zgebal(JOB, N, A, LDA, ILO, IHI, SCALE, INFO)
ZGEBAL
Definition: zgebal.f:162
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine zhseqr(JOB, COMPZ, N, ILO, IHI, H, LDH, W, Z, LDZ, WORK, LWORK, INFO)
ZHSEQR
Definition: zhseqr.f:301
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
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:145
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54
subroutine ztrevc3(SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, LWORK, RWORK, LRWORK, INFO)
ZTREVC3
Definition: ztrevc3.f:248

Here is the call graph for this function:

Here is the caller graph for this function: