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

◆ sgeev()

subroutine sgeev ( character jobvl,
character jobvr,
integer n,
real, dimension( lda, * ) a,
integer lda,
real, dimension( * ) wr,
real, dimension( * ) wi,
real, dimension( ldvl, * ) vl,
integer ldvl,
real, dimension( ldvr, * ) vr,
integer ldvr,
real, dimension( * ) work,
integer lwork,
integer info )

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

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

Purpose:
!>
!> SGEEV computes for an N-by-N real 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 A 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 REAL 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]WR
!>          WR is REAL array, dimension (N)
!> 
[out]WI
!>          WI is REAL array, dimension (N)
!>          WR and WI contain the real and imaginary parts,
!>          respectively, of the computed eigenvalues.  Complex
!>          conjugate pairs of eigenvalues appear consecutively
!>          with the eigenvalue having the positive imaginary part
!>          first.
!> 
[out]VL
!>          VL is REAL 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.
!>          If the j-th eigenvalue is real, then u(j) = VL(:,j),
!>          the j-th column of VL.
!>          If the j-th and (j+1)-st eigenvalues form a complex
!>          conjugate pair, then u(j) = VL(:,j) + i*VL(:,j+1) and
!>          u(j+1) = VL(:,j) - i*VL(:,j+1).
!> 
[in]LDVL
!>          LDVL is INTEGER
!>          The leading dimension of the array VL.  LDVL >= 1; if
!>          JOBVL = 'V', LDVL >= N.
!> 
[out]VR
!>          VR is REAL 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.
!>          If the j-th eigenvalue is real, then v(j) = VR(:,j),
!>          the j-th column of VR.
!>          If the j-th and (j+1)-st eigenvalues form a complex
!>          conjugate pair, then v(j) = VR(:,j) + i*VR(:,j+1) and
!>          v(j+1) = VR(:,j) - i*VR(:,j+1).
!> 
[in]LDVR
!>          LDVR is INTEGER
!>          The leading dimension of the array VR.  LDVR >= 1; if
!>          JOBVR = 'V', LDVR >= N.
!> 
[out]WORK
!>          WORK is REAL 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,3*N), and
!>          if JOBVL = 'V' or JOBVR = 'V', LWORK >= 4*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]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 WR and WI contain eigenvalues which
!>                have converged.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 188 of file sgeev.f.

191 implicit none
192*
193* -- LAPACK driver routine --
194* -- LAPACK is a software package provided by Univ. of Tennessee, --
195* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
196*
197* .. Scalar Arguments ..
198 CHARACTER JOBVL, JOBVR
199 INTEGER INFO, LDA, LDVL, LDVR, LWORK, N
200* ..
201* .. Array Arguments ..
202 REAL A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ),
203 $ WI( * ), WORK( * ), WR( * )
204* ..
205*
206* =====================================================================
207*
208* .. Parameters ..
209 REAL ZERO, ONE
210 parameter( zero = 0.0e0, one = 1.0e0 )
211* ..
212* .. Local Scalars ..
213 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR
214 CHARACTER SIDE
215 INTEGER HSWORK, I, IBAL, IERR, IHI, ILO, ITAU, IWRK, K,
216 $ LWORK_TREVC, MAXWRK, MINWRK, NOUT
217 REAL ANRM, BIGNUM, CS, CSCALE, EPS, R, SCL, SMLNUM,
218 $ SN
219* ..
220* .. Local Arrays ..
221 LOGICAL SELECT( 1 )
222 REAL DUM( 1 )
223* ..
224* .. External Subroutines ..
225 EXTERNAL sgebak, sgebal, sgehrd,
226 $ shseqr, slacpy, slartg,
227 $ slascl, sorghr, srot,
229* ..
230* .. External Functions ..
231 LOGICAL LSAME
232 INTEGER ISAMAX, ILAENV
233 REAL SLAMCH, SLANGE, SLAPY2, SNRM2,
234 $ SROUNDUP_LWORK
235 EXTERNAL lsame, isamax, ilaenv,
236 $ slamch, slange, slapy2,
238* ..
239* .. Intrinsic Functions ..
240 INTRINSIC max, sqrt
241* ..
242* .. Executable Statements ..
243*
244* Test the input arguments
245*
246 info = 0
247 lquery = ( lwork.EQ.-1 )
248 wantvl = lsame( jobvl, 'V' )
249 wantvr = lsame( jobvr, 'V' )
250 IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
251 info = -1
252 ELSE IF( ( .NOT.wantvr ) .AND.
253 $ ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
254 info = -2
255 ELSE IF( n.LT.0 ) THEN
256 info = -3
257 ELSE IF( lda.LT.max( 1, n ) ) THEN
258 info = -5
259 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
260 info = -9
261 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
262 info = -11
263 END IF
264*
265* Compute workspace
266* (Note: Comments in the code beginning "Workspace:" describe the
267* minimal amount of workspace needed at that point in the code,
268* as well as the preferred amount for good performance.
269* NB refers to the optimal block size for the immediately
270* following subroutine, as returned by ILAENV.
271* HSWORK refers to the workspace preferred by SHSEQR, as
272* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
273* the worst case.)
274*
275 IF( info.EQ.0 ) THEN
276 IF( n.EQ.0 ) THEN
277 minwrk = 1
278 maxwrk = 1
279 ELSE
280 maxwrk = 2*n + n*ilaenv( 1, 'SGEHRD', ' ', n, 1, n, 0 )
281 IF( wantvl ) THEN
282 minwrk = 4*n
283 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
284 $ 'SORGHR', ' ', n, 1, n, -1 ) )
285 CALL shseqr( 'S', 'V', n, 1, n, a, lda, wr, wi, vl,
286 $ ldvl,
287 $ work, -1, info )
288 hswork = int( work(1) )
289 maxwrk = max( maxwrk, n + 1, n + hswork )
290 CALL strevc3( 'L', 'B', SELECT, n, a, lda,
291 $ vl, ldvl, vr, ldvr, n, nout,
292 $ work, -1, ierr )
293 lwork_trevc = int( work(1) )
294 maxwrk = max( maxwrk, n + lwork_trevc )
295 maxwrk = max( maxwrk, 4*n )
296 ELSE IF( wantvr ) THEN
297 minwrk = 4*n
298 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
299 $ 'SORGHR', ' ', n, 1, n, -1 ) )
300 CALL shseqr( 'S', 'V', n, 1, n, a, lda, wr, wi, vr,
301 $ ldvr,
302 $ work, -1, info )
303 hswork = int( work(1) )
304 maxwrk = max( maxwrk, n + 1, n + hswork )
305 CALL strevc3( 'R', 'B', SELECT, n, a, lda,
306 $ vl, ldvl, vr, ldvr, n, nout,
307 $ work, -1, ierr )
308 lwork_trevc = int( work(1) )
309 maxwrk = max( maxwrk, n + lwork_trevc )
310 maxwrk = max( maxwrk, 4*n )
311 ELSE
312 minwrk = 3*n
313 CALL shseqr( 'E', 'N', n, 1, n, a, lda, wr, wi, vr,
314 $ ldvr,
315 $ work, -1, info )
316 hswork = int( work(1) )
317 maxwrk = max( maxwrk, n + 1, n + hswork )
318 END IF
319 maxwrk = max( maxwrk, minwrk )
320 END IF
321 work( 1 ) = sroundup_lwork(maxwrk)
322*
323 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
324 info = -13
325 END IF
326 END IF
327*
328 IF( info.NE.0 ) THEN
329 CALL xerbla( 'SGEEV ', -info )
330 RETURN
331 ELSE IF( lquery ) THEN
332 RETURN
333 END IF
334*
335* Quick return if possible
336*
337 IF( n.EQ.0 )
338 $ RETURN
339*
340* Get machine constants
341*
342 eps = slamch( 'P' )
343 smlnum = slamch( 'S' )
344 bignum = one / smlnum
345 smlnum = sqrt( smlnum ) / eps
346 bignum = one / smlnum
347*
348* Scale A if max element outside range [SMLNUM,BIGNUM]
349*
350 anrm = slange( 'M', n, n, a, lda, dum )
351 scalea = .false.
352 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
353 scalea = .true.
354 cscale = smlnum
355 ELSE IF( anrm.GT.bignum ) THEN
356 scalea = .true.
357 cscale = bignum
358 END IF
359 IF( scalea )
360 $ CALL slascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
361*
362* Balance the matrix
363* (Workspace: need N)
364*
365 ibal = 1
366 CALL sgebal( 'B', n, a, lda, ilo, ihi, work( ibal ), ierr )
367*
368* Reduce to upper Hessenberg form
369* (Workspace: need 3*N, prefer 2*N+N*NB)
370*
371 itau = ibal + n
372 iwrk = itau + n
373 CALL sgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
374 $ lwork-iwrk+1, ierr )
375*
376 IF( wantvl ) THEN
377*
378* Want left eigenvectors
379* Copy Householder vectors to VL
380*
381 side = 'L'
382 CALL slacpy( 'L', n, n, a, lda, vl, ldvl )
383*
384* Generate orthogonal matrix in VL
385* (Workspace: need 3*N-1, prefer 2*N+(N-1)*NB)
386*
387 CALL sorghr( n, ilo, ihi, vl, ldvl, work( itau ),
388 $ work( iwrk ),
389 $ lwork-iwrk+1, ierr )
390*
391* Perform QR iteration, accumulating Schur vectors in VL
392* (Workspace: need N+1, prefer N+HSWORK (see comments) )
393*
394 iwrk = itau
395 CALL shseqr( 'S', 'V', n, ilo, ihi, a, lda, wr, wi, vl,
396 $ ldvl,
397 $ work( iwrk ), lwork-iwrk+1, info )
398*
399 IF( wantvr ) THEN
400*
401* Want left and right eigenvectors
402* Copy Schur vectors to VR
403*
404 side = 'B'
405 CALL slacpy( 'F', n, n, vl, ldvl, vr, ldvr )
406 END IF
407*
408 ELSE IF( wantvr ) THEN
409*
410* Want right eigenvectors
411* Copy Householder vectors to VR
412*
413 side = 'R'
414 CALL slacpy( 'L', n, n, a, lda, vr, ldvr )
415*
416* Generate orthogonal matrix in VR
417* (Workspace: need 3*N-1, prefer 2*N+(N-1)*NB)
418*
419 CALL sorghr( n, ilo, ihi, vr, ldvr, work( itau ),
420 $ work( iwrk ),
421 $ lwork-iwrk+1, ierr )
422*
423* Perform QR iteration, accumulating Schur vectors in VR
424* (Workspace: need N+1, prefer N+HSWORK (see comments) )
425*
426 iwrk = itau
427 CALL shseqr( 'S', 'V', n, ilo, ihi, a, lda, wr, wi, vr,
428 $ ldvr,
429 $ work( iwrk ), lwork-iwrk+1, info )
430*
431 ELSE
432*
433* Compute eigenvalues only
434* (Workspace: need N+1, prefer N+HSWORK (see comments) )
435*
436 iwrk = itau
437 CALL shseqr( 'E', 'N', n, ilo, ihi, a, lda, wr, wi, vr,
438 $ ldvr,
439 $ work( iwrk ), lwork-iwrk+1, info )
440 END IF
441*
442* If INFO .NE. 0 from SHSEQR, then quit
443*
444 IF( info.NE.0 )
445 $ GO TO 50
446*
447 IF( wantvl .OR. wantvr ) THEN
448*
449* Compute left and/or right eigenvectors
450* (Workspace: need 4*N, prefer N + N + 2*N*NB)
451*
452 CALL strevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr,
453 $ ldvr,
454 $ n, nout, work( iwrk ), lwork-iwrk+1, ierr )
455 END IF
456*
457 IF( wantvl ) THEN
458*
459* Undo balancing of left eigenvectors
460* (Workspace: need N)
461*
462 CALL sgebak( 'B', 'L', n, ilo, ihi, work( ibal ), n, vl,
463 $ ldvl,
464 $ ierr )
465*
466* Normalize left eigenvectors and make largest component real
467*
468 DO 20 i = 1, n
469 IF( wi( i ).EQ.zero ) THEN
470 scl = one / snrm2( n, vl( 1, i ), 1 )
471 CALL sscal( n, scl, vl( 1, i ), 1 )
472 ELSE IF( wi( i ).GT.zero ) THEN
473 scl = one / slapy2( snrm2( n, vl( 1, i ), 1 ),
474 $ snrm2( n, vl( 1, i+1 ), 1 ) )
475 CALL sscal( n, scl, vl( 1, i ), 1 )
476 CALL sscal( n, scl, vl( 1, i+1 ), 1 )
477 DO 10 k = 1, n
478 work( iwrk+k-1 ) = vl( k, i )**2 + vl( k, i+1 )**2
479 10 CONTINUE
480 k = isamax( n, work( iwrk ), 1 )
481 CALL slartg( vl( k, i ), vl( k, i+1 ), cs, sn, r )
482 CALL srot( n, vl( 1, i ), 1, vl( 1, i+1 ), 1, cs, sn )
483 vl( k, i+1 ) = zero
484 END IF
485 20 CONTINUE
486 END IF
487*
488 IF( wantvr ) THEN
489*
490* Undo balancing of right eigenvectors
491* (Workspace: need N)
492*
493 CALL sgebak( 'B', 'R', n, ilo, ihi, work( ibal ), n, vr,
494 $ ldvr,
495 $ ierr )
496*
497* Normalize right eigenvectors and make largest component real
498*
499 DO 40 i = 1, n
500 IF( wi( i ).EQ.zero ) THEN
501 scl = one / snrm2( n, vr( 1, i ), 1 )
502 CALL sscal( n, scl, vr( 1, i ), 1 )
503 ELSE IF( wi( i ).GT.zero ) THEN
504 scl = one / slapy2( snrm2( n, vr( 1, i ), 1 ),
505 $ snrm2( n, vr( 1, i+1 ), 1 ) )
506 CALL sscal( n, scl, vr( 1, i ), 1 )
507 CALL sscal( n, scl, vr( 1, i+1 ), 1 )
508 DO 30 k = 1, n
509 work( iwrk+k-1 ) = vr( k, i )**2 + vr( k, i+1 )**2
510 30 CONTINUE
511 k = isamax( n, work( iwrk ), 1 )
512 CALL slartg( vr( k, i ), vr( k, i+1 ), cs, sn, r )
513 CALL srot( n, vr( 1, i ), 1, vr( 1, i+1 ), 1, cs, sn )
514 vr( k, i+1 ) = zero
515 END IF
516 40 CONTINUE
517 END IF
518*
519* Undo scaling if necessary
520*
521 50 CONTINUE
522 IF( scalea ) THEN
523 CALL slascl( 'G', 0, 0, cscale, anrm, n-info, 1,
524 $ wr( info+1 ),
525 $ max( n-info, 1 ), ierr )
526 CALL slascl( 'G', 0, 0, cscale, anrm, n-info, 1,
527 $ wi( info+1 ),
528 $ max( n-info, 1 ), ierr )
529 IF( info.GT.0 ) THEN
530 CALL slascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wr, n,
531 $ ierr )
532 CALL slascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
533 $ ierr )
534 END IF
535 END IF
536*
537 work( 1 ) = sroundup_lwork(maxwrk)
538 RETURN
539*
540* End of SGEEV
541*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
SGEBAK
Definition sgebak.f:128
subroutine sgebal(job, n, a, lda, ilo, ihi, scale, info)
SGEBAL
Definition sgebal.f:161
subroutine sgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
SGEHRD
Definition sgehrd.f:166
subroutine shseqr(job, compz, n, ilo, ihi, h, ldh, wr, wi, z, ldz, work, lwork, info)
SHSEQR
Definition shseqr.f:314
integer function isamax(n, sx, incx)
ISAMAX
Definition isamax.f:71
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slange(norm, m, n, a, lda, work)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition slange.f:112
real function slapy2(x, y)
SLAPY2 returns sqrt(x2+y2).
Definition slapy2.f:61
subroutine slartg(f, g, c, s, r)
SLARTG generates a plane rotation with real cosine and real sine.
Definition slartg.f90:111
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function snrm2(n, x, incx)
SNRM2
Definition snrm2.f90:89
subroutine srot(n, sx, incx, sy, incy, c, s)
SROT
Definition srot.f:92
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine strevc3(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, lwork, info)
STREVC3
Definition strevc3.f:235
subroutine sorghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
SORGHR
Definition sorghr.f:125
Here is the call graph for this function:
Here is the caller graph for this function: