LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sgegs ( character  JOBVSL,
character  JOBVSR,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  ALPHAR,
real, dimension( * )  ALPHAI,
real, dimension( * )  BETA,
real, dimension( ldvsl, * )  VSL,
integer  LDVSL,
real, dimension( ldvsr, * )  VSR,
integer  LDVSR,
real, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

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

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

Purpose:
 This routine is deprecated and has been replaced by routine SGGES.

 SGEGS computes the eigenvalues, real Schur form, and, optionally,
 left and or/right Schur vectors of a real matrix pair (A,B).
 Given two square matrices A and B, the generalized real Schur
 factorization has the form
 
   A = Q*S*Z**T,  B = Q*T*Z**T

 where Q and Z are orthogonal matrices, T is upper triangular, and S
 is an upper quasi-triangular matrix with 1-by-1 and 2-by-2 diagonal
 blocks, the 2-by-2 blocks corresponding to complex conjugate pairs
 of eigenvalues of (A,B).  The columns of Q are the left Schur vectors
 and the columns of Z are the right Schur vectors.
 
 If only the eigenvalues of (A,B) are needed, the driver routine
 SGEGV should be used instead.  See SGEGV for a description of the
 eigenvalues of the generalized nonsymmetric eigenvalue problem
 (GNEP).
Parameters
[in]JOBVSL
          JOBVSL is CHARACTER*1
          = 'N':  do not compute the left Schur vectors;
          = 'V':  compute the left Schur vectors (returned in VSL).
[in]JOBVSR
          JOBVSR is CHARACTER*1
          = 'N':  do not compute the right Schur vectors;
          = 'V':  compute the right Schur vectors (returned in VSR).
[in]N
          N is INTEGER
          The order of the matrices A, B, VSL, and VSR.  N >= 0.
[in,out]A
          A is REAL array, dimension (LDA, N)
          On entry, the matrix A.
          On exit, the upper quasi-triangular matrix S from the
          generalized real Schur factorization.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension (LDB, N)
          On entry, the matrix B.
          On exit, the upper triangular matrix T from the generalized
          real Schur factorization.
[in]LDB
          LDB is INTEGER
          The leading dimension of B.  LDB >= max(1,N).
[out]ALPHAR
          ALPHAR is REAL array, dimension (N)
          The real parts of each scalar alpha defining an eigenvalue
          of GNEP.
[out]ALPHAI
          ALPHAI is REAL array, dimension (N)
          The imaginary parts of each scalar alpha defining an
          eigenvalue of GNEP.  If ALPHAI(j) is zero, then the j-th
          eigenvalue is real; if positive, then the j-th and (j+1)-st
          eigenvalues are a complex conjugate pair, with
          ALPHAI(j+1) = -ALPHAI(j).
[out]BETA
          BETA is REAL array, dimension (N)
          The scalars beta that define the eigenvalues of GNEP.
          Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and
          beta = BETA(j) represent the j-th eigenvalue of the matrix
          pair (A,B), in one of the forms lambda = alpha/beta or
          mu = beta/alpha.  Since either lambda or mu may overflow,
          they should not, in general, be computed.
[out]VSL
          VSL is REAL array, dimension (LDVSL,N)
          If JOBVSL = 'V', the matrix of left Schur vectors Q.
          Not referenced if JOBVSL = 'N'.
[in]LDVSL
          LDVSL is INTEGER
          The leading dimension of the matrix VSL. LDVSL >=1, and
          if JOBVSL = 'V', LDVSL >= N.
[out]VSR
          VSR is REAL array, dimension (LDVSR,N)
          If JOBVSR = 'V', the matrix of right Schur vectors Z.
          Not referenced if JOBVSR = 'N'.
[in]LDVSR
          LDVSR is INTEGER
          The leading dimension of the matrix VSR. LDVSR >= 1, and
          if JOBVSR = 'V', LDVSR >= 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,4*N).
          For good performance, LWORK must generally be larger.
          To compute the optimal value of LWORK, call ILAENV to get
          blocksizes (for SGEQRF, SORMQR, and SORGQR.)  Then compute:
          NB  -- MAX of the blocksizes for SGEQRF, SORMQR, and SORGQR
          The optimal LWORK is  2*N + N*(NB+1).

          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.
          = 1,...,N:
                The QZ iteration failed.  (A,B) are not in Schur
                form, but ALPHAR(j), ALPHAI(j), and BETA(j) should
                be correct for j=INFO+1,...,N.
          > N:  errors that usually indicate LAPACK problems:
                =N+1: error return from SGGBAL
                =N+2: error return from SGEQRF
                =N+3: error return from SORMQR
                =N+4: error return from SORGQR
                =N+5: error return from SGGHRD
                =N+6: error return from SHGEQZ (other than failed
                                                iteration)
                =N+7: error return from SGGBAK (computing VSL)
                =N+8: error return from SGGBAK (computing VSR)
                =N+9: error return from SLASCL (various places)
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 229 of file sgegs.f.

229 *
230 * -- LAPACK driver routine (version 3.4.0) --
231 * -- LAPACK is a software package provided by Univ. of Tennessee, --
232 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
233 * November 2011
234 *
235 * .. Scalar Arguments ..
236  CHARACTER jobvsl, jobvsr
237  INTEGER info, lda, ldb, ldvsl, ldvsr, lwork, n
238 * ..
239 * .. Array Arguments ..
240  REAL a( lda, * ), alphai( * ), alphar( * ),
241  $ b( ldb, * ), beta( * ), vsl( ldvsl, * ),
242  $ vsr( ldvsr, * ), work( * )
243 * ..
244 *
245 * =====================================================================
246 *
247 * .. Parameters ..
248  REAL zero, one
249  parameter ( zero = 0.0e0, one = 1.0e0 )
250 * ..
251 * .. Local Scalars ..
252  LOGICAL ilascl, ilbscl, ilvsl, ilvsr, lquery
253  INTEGER icols, ihi, iinfo, ijobvl, ijobvr, ileft,
254  $ ilo, iright, irows, itau, iwork, lopt, lwkmin,
255  $ lwkopt, nb, nb1, nb2, nb3
256  REAL anrm, anrmto, bignum, bnrm, bnrmto, eps,
257  $ safmin, smlnum
258 * ..
259 * .. External Subroutines ..
260  EXTERNAL sgeqrf, sggbak, sggbal, sgghrd, shgeqz, slacpy,
262 * ..
263 * .. External Functions ..
264  LOGICAL lsame
265  INTEGER ilaenv
266  REAL slamch, slange
267  EXTERNAL ilaenv, lsame, slamch, slange
268 * ..
269 * .. Intrinsic Functions ..
270  INTRINSIC int, max
271 * ..
272 * .. Executable Statements ..
273 *
274 * Decode the input arguments
275 *
276  IF( lsame( jobvsl, 'N' ) ) THEN
277  ijobvl = 1
278  ilvsl = .false.
279  ELSE IF( lsame( jobvsl, 'V' ) ) THEN
280  ijobvl = 2
281  ilvsl = .true.
282  ELSE
283  ijobvl = -1
284  ilvsl = .false.
285  END IF
286 *
287  IF( lsame( jobvsr, 'N' ) ) THEN
288  ijobvr = 1
289  ilvsr = .false.
290  ELSE IF( lsame( jobvsr, 'V' ) ) THEN
291  ijobvr = 2
292  ilvsr = .true.
293  ELSE
294  ijobvr = -1
295  ilvsr = .false.
296  END IF
297 *
298 * Test the input arguments
299 *
300  lwkmin = max( 4*n, 1 )
301  lwkopt = lwkmin
302  work( 1 ) = lwkopt
303  lquery = ( lwork.EQ.-1 )
304  info = 0
305  IF( ijobvl.LE.0 ) THEN
306  info = -1
307  ELSE IF( ijobvr.LE.0 ) THEN
308  info = -2
309  ELSE IF( n.LT.0 ) THEN
310  info = -3
311  ELSE IF( lda.LT.max( 1, n ) ) THEN
312  info = -5
313  ELSE IF( ldb.LT.max( 1, n ) ) THEN
314  info = -7
315  ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
316  info = -12
317  ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
318  info = -14
319  ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
320  info = -16
321  END IF
322 *
323  IF( info.EQ.0 ) THEN
324  nb1 = ilaenv( 1, 'SGEQRF', ' ', n, n, -1, -1 )
325  nb2 = ilaenv( 1, 'SORMQR', ' ', n, n, n, -1 )
326  nb3 = ilaenv( 1, 'SORGQR', ' ', n, n, n, -1 )
327  nb = max( nb1, nb2, nb3 )
328  lopt = 2*n+n*(nb+1)
329  work( 1 ) = lopt
330  END IF
331 *
332  IF( info.NE.0 ) THEN
333  CALL xerbla( 'SGEGS ', -info )
334  RETURN
335  ELSE IF( lquery ) THEN
336  RETURN
337  END IF
338 *
339 * Quick return if possible
340 *
341  IF( n.EQ.0 )
342  $ RETURN
343 *
344 * Get machine constants
345 *
346  eps = slamch( 'E' )*slamch( 'B' )
347  safmin = slamch( 'S' )
348  smlnum = n*safmin / eps
349  bignum = one / smlnum
350 *
351 * Scale A if max element outside range [SMLNUM,BIGNUM]
352 *
353  anrm = slange( 'M', n, n, a, lda, work )
354  ilascl = .false.
355  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
356  anrmto = smlnum
357  ilascl = .true.
358  ELSE IF( anrm.GT.bignum ) THEN
359  anrmto = bignum
360  ilascl = .true.
361  END IF
362 *
363  IF( ilascl ) THEN
364  CALL slascl( 'G', -1, -1, anrm, anrmto, n, n, a, lda, iinfo )
365  IF( iinfo.NE.0 ) THEN
366  info = n + 9
367  RETURN
368  END IF
369  END IF
370 *
371 * Scale B if max element outside range [SMLNUM,BIGNUM]
372 *
373  bnrm = slange( 'M', n, n, b, ldb, work )
374  ilbscl = .false.
375  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
376  bnrmto = smlnum
377  ilbscl = .true.
378  ELSE IF( bnrm.GT.bignum ) THEN
379  bnrmto = bignum
380  ilbscl = .true.
381  END IF
382 *
383  IF( ilbscl ) THEN
384  CALL slascl( 'G', -1, -1, bnrm, bnrmto, n, n, b, ldb, iinfo )
385  IF( iinfo.NE.0 ) THEN
386  info = n + 9
387  RETURN
388  END IF
389  END IF
390 *
391 * Permute the matrix to make it more nearly triangular
392 * Workspace layout: (2*N words -- "work..." not actually used)
393 * left_permutation, right_permutation, work...
394 *
395  ileft = 1
396  iright = n + 1
397  iwork = iright + n
398  CALL sggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
399  $ work( iright ), work( iwork ), iinfo )
400  IF( iinfo.NE.0 ) THEN
401  info = n + 1
402  GO TO 10
403  END IF
404 *
405 * Reduce B to triangular form, and initialize VSL and/or VSR
406 * Workspace layout: ("work..." must have at least N words)
407 * left_permutation, right_permutation, tau, work...
408 *
409  irows = ihi + 1 - ilo
410  icols = n + 1 - ilo
411  itau = iwork
412  iwork = itau + irows
413  CALL sgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
414  $ work( iwork ), lwork+1-iwork, iinfo )
415  IF( iinfo.GE.0 )
416  $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
417  IF( iinfo.NE.0 ) THEN
418  info = n + 2
419  GO TO 10
420  END IF
421 *
422  CALL sormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
423  $ work( itau ), a( ilo, ilo ), lda, work( iwork ),
424  $ lwork+1-iwork, iinfo )
425  IF( iinfo.GE.0 )
426  $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
427  IF( iinfo.NE.0 ) THEN
428  info = n + 3
429  GO TO 10
430  END IF
431 *
432  IF( ilvsl ) THEN
433  CALL slaset( 'Full', n, n, zero, one, vsl, ldvsl )
434  CALL slacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
435  $ vsl( ilo+1, ilo ), ldvsl )
436  CALL sorgqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
437  $ work( itau ), work( iwork ), lwork+1-iwork,
438  $ iinfo )
439  IF( iinfo.GE.0 )
440  $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
441  IF( iinfo.NE.0 ) THEN
442  info = n + 4
443  GO TO 10
444  END IF
445  END IF
446 *
447  IF( ilvsr )
448  $ CALL slaset( 'Full', n, n, zero, one, vsr, ldvsr )
449 *
450 * Reduce to generalized Hessenberg form
451 *
452  CALL sgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
453  $ ldvsl, vsr, ldvsr, iinfo )
454  IF( iinfo.NE.0 ) THEN
455  info = n + 5
456  GO TO 10
457  END IF
458 *
459 * Perform QZ algorithm, computing Schur vectors if desired
460 * Workspace layout: ("work..." must have at least 1 word)
461 * left_permutation, right_permutation, work...
462 *
463  iwork = itau
464  CALL shgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
465  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
466  $ work( iwork ), lwork+1-iwork, iinfo )
467  IF( iinfo.GE.0 )
468  $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
469  IF( iinfo.NE.0 ) THEN
470  IF( iinfo.GT.0 .AND. iinfo.LE.n ) THEN
471  info = iinfo
472  ELSE IF( iinfo.GT.n .AND. iinfo.LE.2*n ) THEN
473  info = iinfo - n
474  ELSE
475  info = n + 6
476  END IF
477  GO TO 10
478  END IF
479 *
480 * Apply permutation to VSL and VSR
481 *
482  IF( ilvsl ) THEN
483  CALL sggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
484  $ work( iright ), n, vsl, ldvsl, iinfo )
485  IF( iinfo.NE.0 ) THEN
486  info = n + 7
487  GO TO 10
488  END IF
489  END IF
490  IF( ilvsr ) THEN
491  CALL sggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
492  $ work( iright ), n, vsr, ldvsr, iinfo )
493  IF( iinfo.NE.0 ) THEN
494  info = n + 8
495  GO TO 10
496  END IF
497  END IF
498 *
499 * Undo scaling
500 *
501  IF( ilascl ) THEN
502  CALL slascl( 'H', -1, -1, anrmto, anrm, n, n, a, lda, iinfo )
503  IF( iinfo.NE.0 ) THEN
504  info = n + 9
505  RETURN
506  END IF
507  CALL slascl( 'G', -1, -1, anrmto, anrm, n, 1, alphar, n,
508  $ iinfo )
509  IF( iinfo.NE.0 ) THEN
510  info = n + 9
511  RETURN
512  END IF
513  CALL slascl( 'G', -1, -1, anrmto, anrm, n, 1, alphai, n,
514  $ iinfo )
515  IF( iinfo.NE.0 ) THEN
516  info = n + 9
517  RETURN
518  END IF
519  END IF
520 *
521  IF( ilbscl ) THEN
522  CALL slascl( 'U', -1, -1, bnrmto, bnrm, n, n, b, ldb, iinfo )
523  IF( iinfo.NE.0 ) THEN
524  info = n + 9
525  RETURN
526  END IF
527  CALL slascl( 'G', -1, -1, bnrmto, bnrm, n, 1, beta, n, iinfo )
528  IF( iinfo.NE.0 ) THEN
529  info = n + 9
530  RETURN
531  END IF
532  END IF
533 *
534  10 CONTINUE
535  work( 1 ) = lwkopt
536 *
537  RETURN
538 *
539 * End of SGEGS
540 *
subroutine sggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
SGGBAL
Definition: sggbal.f:179
subroutine sormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMQR
Definition: sormqr.f:170
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:145
subroutine sgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
SGEQRF
Definition: sgeqrf.f:138
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine shgeqz(JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, INFO)
SHGEQZ
Definition: shgeqz.f:306
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
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:116
subroutine sggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
SGGBAK
Definition: sggbak.f:149
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine sgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
SGGHRD
Definition: sgghrd.f:209
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine sorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
SORGQR
Definition: sorgqr.f:130
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function: