LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ sggesx()

subroutine sggesx ( character  JOBVSL,
character  JOBVSR,
character  SORT,
external  SELCTG,
character  SENSE,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
integer  SDIM,
real, dimension( * )  ALPHAR,
real, dimension( * )  ALPHAI,
real, dimension( * )  BETA,
real, dimension( ldvsl, * )  VSL,
integer  LDVSL,
real, dimension( ldvsr, * )  VSR,
integer  LDVSR,
real, dimension( 2 )  RCONDE,
real, dimension( 2 )  RCONDV,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
logical, dimension( * )  BWORK,
integer  INFO 
)

SGGESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices

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

Purpose:
 SGGESX computes for a pair of N-by-N real nonsymmetric matrices
 (A,B), the generalized eigenvalues, the real Schur form (S,T), and,
 optionally, the left and/or right matrices of Schur vectors (VSL and
 VSR).  This gives the generalized Schur factorization

      (A,B) = ( (VSL) S (VSR)**T, (VSL) T (VSR)**T )

 Optionally, it also orders the eigenvalues so that a selected cluster
 of eigenvalues appears in the leading diagonal blocks of the upper
 quasi-triangular matrix S and the upper triangular matrix T; computes
 a reciprocal condition number for the average of the selected
 eigenvalues (RCONDE); and computes a reciprocal condition number for
 the right and left deflating subspaces corresponding to the selected
 eigenvalues (RCONDV). The leading columns of VSL and VSR then form
 an orthonormal basis for the corresponding left and right eigenspaces
 (deflating subspaces).

 A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
 or a ratio alpha/beta = w, such that  A - w*B is singular.  It is
 usually represented as the pair (alpha,beta), as there is a
 reasonable interpretation for beta=0 or for both being zero.

 A pair of matrices (S,T) is in generalized real Schur form if T is
 upper triangular with non-negative diagonal and S is block upper
 triangular with 1-by-1 and 2-by-2 blocks.  1-by-1 blocks correspond
 to real generalized eigenvalues, while 2-by-2 blocks of S will be
 "standardized" by making the corresponding elements of T have the
 form:
         [  a  0  ]
         [  0  b  ]

 and the pair of corresponding 2-by-2 blocks in S and T will have a
 complex conjugate pair of generalized eigenvalues.
Parameters
[in]JOBVSL
          JOBVSL is CHARACTER*1
          = 'N':  do not compute the left Schur vectors;
          = 'V':  compute the left Schur vectors.
[in]JOBVSR
          JOBVSR is CHARACTER*1
          = 'N':  do not compute the right Schur vectors;
          = 'V':  compute the right Schur vectors.
[in]SORT
          SORT is CHARACTER*1
          Specifies whether or not to order the eigenvalues on the
          diagonal of the generalized Schur form.
          = 'N':  Eigenvalues are not ordered;
          = 'S':  Eigenvalues are ordered (see SELCTG).
[in]SELCTG
          SELCTG is a LOGICAL FUNCTION of three REAL arguments
          SELCTG must be declared EXTERNAL in the calling subroutine.
          If SORT = 'N', SELCTG is not referenced.
          If SORT = 'S', SELCTG is used to select eigenvalues to sort
          to the top left of the Schur form.
          An eigenvalue (ALPHAR(j)+ALPHAI(j))/BETA(j) is selected if
          SELCTG(ALPHAR(j),ALPHAI(j),BETA(j)) is true; i.e. if either
          one of a complex conjugate pair of eigenvalues is selected,
          then both complex eigenvalues are selected.
          Note that a selected complex eigenvalue may no longer satisfy
          SELCTG(ALPHAR(j),ALPHAI(j),BETA(j)) = .TRUE. after ordering,
          since ordering may change the value of complex eigenvalues
          (especially if the eigenvalue is ill-conditioned), in this
          case INFO is set to N+3.
[in]SENSE
          SENSE is CHARACTER*1
          Determines which reciprocal condition numbers are computed.
          = 'N':  None are computed;
          = 'E':  Computed for average of selected eigenvalues only;
          = 'V':  Computed for selected deflating subspaces only;
          = 'B':  Computed for both.
          If SENSE = 'E', 'V', or 'B', SORT must equal 'S'.
[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 first of the pair of matrices.
          On exit, A has been overwritten by its generalized Schur
          form S.
[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 second of the pair of matrices.
          On exit, B has been overwritten by its generalized Schur
          form T.
[in]LDB
          LDB is INTEGER
          The leading dimension of B.  LDB >= max(1,N).
[out]SDIM
          SDIM is INTEGER
          If SORT = 'N', SDIM = 0.
          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
          for which SELCTG is true.  (Complex conjugate pairs for which
          SELCTG is true for either eigenvalue count as 2.)
[out]ALPHAR
          ALPHAR is REAL array, dimension (N)
[out]ALPHAI
          ALPHAI is REAL array, dimension (N)
[out]BETA
          BETA is REAL array, dimension (N)
          On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
          be the generalized eigenvalues.  ALPHAR(j) + ALPHAI(j)*i
          and BETA(j),j=1,...,N  are the diagonals of the complex Schur
          form (S,T) that would result if the 2-by-2 diagonal blocks of
          the real Schur form of (A,B) were further reduced to
          triangular form using 2-by-2 complex unitary transformations.
          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) negative.

          Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j)
          may easily over- or underflow, and BETA(j) may even be zero.
          Thus, the user should avoid naively computing the ratio.
          However, ALPHAR and ALPHAI will be always less than and
          usually comparable with norm(A) in magnitude, and BETA always
          less than and usually comparable with norm(B).
[out]VSL
          VSL is REAL array, dimension (LDVSL,N)
          If JOBVSL = 'V', VSL will contain the left Schur vectors.
          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', VSR will contain the right Schur vectors.
          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]RCONDE
          RCONDE is REAL array, dimension ( 2 )
          If SENSE = 'E' or 'B', RCONDE(1) and RCONDE(2) contain the
          reciprocal condition numbers for the average of the selected
          eigenvalues.
          Not referenced if SENSE = 'N' or 'V'.
[out]RCONDV
          RCONDV is REAL array, dimension ( 2 )
          If SENSE = 'V' or 'B', RCONDV(1) and RCONDV(2) contain the
          reciprocal condition numbers for the selected deflating
          subspaces.
          Not referenced if SENSE = 'N' or 'E'.
[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.
          If N = 0, LWORK >= 1, else if SENSE = 'E', 'V', or 'B',
          LWORK >= max( 8*N, 6*N+16, 2*SDIM*(N-SDIM) ), else
          LWORK >= max( 8*N, 6*N+16 ).
          Note that 2*SDIM*(N-SDIM) <= N*N/2.
          Note also that an error is only returned if
          LWORK < max( 8*N, 6*N+16), but if SENSE = 'E' or 'V' or 'B'
          this may not be large enough.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the bound on the optimal size of the WORK
          array and the minimum size of the IWORK array, returns these
          values as the first entries of the WORK and IWORK arrays, and
          no error message related to LWORK or LIWORK is issued by
          XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If SENSE = 'N' or N = 0, LIWORK >= 1, otherwise
          LIWORK >= N+6.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the bound on the optimal size of the
          WORK array and the minimum size of the IWORK array, returns
          these values as the first entries of the WORK and IWORK
          arrays, and no error message related to LWORK or LIWORK is
          issued by XERBLA.
[out]BWORK
          BWORK is LOGICAL array, dimension (N)
          Not referenced if SORT = 'N'.
[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:  =N+1: other than QZ iteration failed in SHGEQZ
                =N+2: after reordering, roundoff changed values of
                      some complex eigenvalues so that leading
                      eigenvalues in the Generalized Schur form no
                      longer satisfy SELCTG=.TRUE.  This could also
                      be caused due to scaling.
                =N+3: reordering failed in STGSEN.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  An approximate (asymptotic) bound on the average absolute error of
  the selected eigenvalues is

       EPS * norm((A, B)) / RCONDE( 1 ).

  An approximate (asymptotic) bound on the maximum angular error in
  the computed deflating subspaces is

       EPS * norm((A, B)) / RCONDV( 2 ).

  See LAPACK User's Guide, section 4.11 for more information.

Definition at line 361 of file sggesx.f.

365 *
366 * -- LAPACK driver routine --
367 * -- LAPACK is a software package provided by Univ. of Tennessee, --
368 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
369 *
370 * .. Scalar Arguments ..
371  CHARACTER JOBVSL, JOBVSR, SENSE, SORT
372  INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LIWORK, LWORK, N,
373  $ SDIM
374 * ..
375 * .. Array Arguments ..
376  LOGICAL BWORK( * )
377  INTEGER IWORK( * )
378  REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
379  $ B( LDB, * ), BETA( * ), RCONDE( 2 ),
380  $ RCONDV( 2 ), VSL( LDVSL, * ), VSR( LDVSR, * ),
381  $ WORK( * )
382 * ..
383 * .. Function Arguments ..
384  LOGICAL SELCTG
385  EXTERNAL selctg
386 * ..
387 *
388 * =====================================================================
389 *
390 * .. Parameters ..
391  REAL ZERO, ONE
392  parameter( zero = 0.0e+0, one = 1.0e+0 )
393 * ..
394 * .. Local Scalars ..
395  LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
396  $ LQUERY, LST2SL, WANTSB, WANTSE, WANTSN, WANTST,
397  $ WANTSV
398  INTEGER I, ICOLS, IERR, IHI, IJOB, IJOBVL, IJOBVR,
399  $ ILEFT, ILO, IP, IRIGHT, IROWS, ITAU, IWRK,
400  $ LIWMIN, LWRK, MAXWRK, MINWRK
401  REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PL,
402  $ PR, SAFMAX, SAFMIN, SMLNUM
403 * ..
404 * .. Local Arrays ..
405  REAL DIF( 2 )
406 * ..
407 * .. External Subroutines ..
408  EXTERNAL sgeqrf, sggbak, sggbal, sgghrd, shgeqz, slabad,
410  $ xerbla
411 * ..
412 * .. External Functions ..
413  LOGICAL LSAME
414  INTEGER ILAENV
415  REAL SLAMCH, SLANGE
416  EXTERNAL lsame, ilaenv, slamch, slange
417 * ..
418 * .. Intrinsic Functions ..
419  INTRINSIC abs, max, sqrt
420 * ..
421 * .. Executable Statements ..
422 *
423 * Decode the input arguments
424 *
425  IF( lsame( jobvsl, 'N' ) ) THEN
426  ijobvl = 1
427  ilvsl = .false.
428  ELSE IF( lsame( jobvsl, 'V' ) ) THEN
429  ijobvl = 2
430  ilvsl = .true.
431  ELSE
432  ijobvl = -1
433  ilvsl = .false.
434  END IF
435 *
436  IF( lsame( jobvsr, 'N' ) ) THEN
437  ijobvr = 1
438  ilvsr = .false.
439  ELSE IF( lsame( jobvsr, 'V' ) ) THEN
440  ijobvr = 2
441  ilvsr = .true.
442  ELSE
443  ijobvr = -1
444  ilvsr = .false.
445  END IF
446 *
447  wantst = lsame( sort, 'S' )
448  wantsn = lsame( sense, 'N' )
449  wantse = lsame( sense, 'E' )
450  wantsv = lsame( sense, 'V' )
451  wantsb = lsame( sense, 'B' )
452  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
453  IF( wantsn ) THEN
454  ijob = 0
455  ELSE IF( wantse ) THEN
456  ijob = 1
457  ELSE IF( wantsv ) THEN
458  ijob = 2
459  ELSE IF( wantsb ) THEN
460  ijob = 4
461  END IF
462 *
463 * Test the input arguments
464 *
465  info = 0
466  IF( ijobvl.LE.0 ) THEN
467  info = -1
468  ELSE IF( ijobvr.LE.0 ) THEN
469  info = -2
470  ELSE IF( ( .NOT.wantst ) .AND. ( .NOT.lsame( sort, 'N' ) ) ) THEN
471  info = -3
472  ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
473  $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
474  info = -5
475  ELSE IF( n.LT.0 ) THEN
476  info = -6
477  ELSE IF( lda.LT.max( 1, n ) ) THEN
478  info = -8
479  ELSE IF( ldb.LT.max( 1, n ) ) THEN
480  info = -10
481  ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
482  info = -16
483  ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
484  info = -18
485  END IF
486 *
487 * Compute workspace
488 * (Note: Comments in the code beginning "Workspace:" describe the
489 * minimal amount of workspace needed at that point in the code,
490 * as well as the preferred amount for good performance.
491 * NB refers to the optimal block size for the immediately
492 * following subroutine, as returned by ILAENV.)
493 *
494  IF( info.EQ.0 ) THEN
495  IF( n.GT.0) THEN
496  minwrk = max( 8*n, 6*n + 16 )
497  maxwrk = minwrk - n +
498  $ n*ilaenv( 1, 'SGEQRF', ' ', n, 1, n, 0 )
499  maxwrk = max( maxwrk, minwrk - n +
500  $ n*ilaenv( 1, 'SORMQR', ' ', n, 1, n, -1 ) )
501  IF( ilvsl ) THEN
502  maxwrk = max( maxwrk, minwrk - n +
503  $ n*ilaenv( 1, 'SORGQR', ' ', n, 1, n, -1 ) )
504  END IF
505  lwrk = maxwrk
506  IF( ijob.GE.1 )
507  $ lwrk = max( lwrk, n*n/2 )
508  ELSE
509  minwrk = 1
510  maxwrk = 1
511  lwrk = 1
512  END IF
513  work( 1 ) = lwrk
514  IF( wantsn .OR. n.EQ.0 ) THEN
515  liwmin = 1
516  ELSE
517  liwmin = n + 6
518  END IF
519  iwork( 1 ) = liwmin
520 *
521  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
522  info = -22
523  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
524  info = -24
525  END IF
526  END IF
527 *
528  IF( info.NE.0 ) THEN
529  CALL xerbla( 'SGGESX', -info )
530  RETURN
531  ELSE IF (lquery) THEN
532  RETURN
533  END IF
534 *
535 * Quick return if possible
536 *
537  IF( n.EQ.0 ) THEN
538  sdim = 0
539  RETURN
540  END IF
541 *
542 * Get machine constants
543 *
544  eps = slamch( 'P' )
545  safmin = slamch( 'S' )
546  safmax = one / safmin
547  CALL slabad( safmin, safmax )
548  smlnum = sqrt( safmin ) / eps
549  bignum = one / smlnum
550 *
551 * Scale A if max element outside range [SMLNUM,BIGNUM]
552 *
553  anrm = slange( 'M', n, n, a, lda, work )
554  ilascl = .false.
555  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
556  anrmto = smlnum
557  ilascl = .true.
558  ELSE IF( anrm.GT.bignum ) THEN
559  anrmto = bignum
560  ilascl = .true.
561  END IF
562  IF( ilascl )
563  $ CALL slascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
564 *
565 * Scale B if max element outside range [SMLNUM,BIGNUM]
566 *
567  bnrm = slange( 'M', n, n, b, ldb, work )
568  ilbscl = .false.
569  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
570  bnrmto = smlnum
571  ilbscl = .true.
572  ELSE IF( bnrm.GT.bignum ) THEN
573  bnrmto = bignum
574  ilbscl = .true.
575  END IF
576  IF( ilbscl )
577  $ CALL slascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
578 *
579 * Permute the matrix to make it more nearly triangular
580 * (Workspace: need 6*N + 2*N for permutation parameters)
581 *
582  ileft = 1
583  iright = n + 1
584  iwrk = iright + n
585  CALL sggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
586  $ work( iright ), work( iwrk ), ierr )
587 *
588 * Reduce B to triangular form (QR decomposition of B)
589 * (Workspace: need N, prefer N*NB)
590 *
591  irows = ihi + 1 - ilo
592  icols = n + 1 - ilo
593  itau = iwrk
594  iwrk = itau + irows
595  CALL sgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
596  $ work( iwrk ), lwork+1-iwrk, ierr )
597 *
598 * Apply the orthogonal transformation to matrix A
599 * (Workspace: need N, prefer N*NB)
600 *
601  CALL sormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
602  $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
603  $ lwork+1-iwrk, ierr )
604 *
605 * Initialize VSL
606 * (Workspace: need N, prefer N*NB)
607 *
608  IF( ilvsl ) THEN
609  CALL slaset( 'Full', n, n, zero, one, vsl, ldvsl )
610  IF( irows.GT.1 ) THEN
611  CALL slacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
612  $ vsl( ilo+1, ilo ), ldvsl )
613  END IF
614  CALL sorgqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
615  $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
616  END IF
617 *
618 * Initialize VSR
619 *
620  IF( ilvsr )
621  $ CALL slaset( 'Full', n, n, zero, one, vsr, ldvsr )
622 *
623 * Reduce to generalized Hessenberg form
624 * (Workspace: none needed)
625 *
626  CALL sgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
627  $ ldvsl, vsr, ldvsr, ierr )
628 *
629  sdim = 0
630 *
631 * Perform QZ algorithm, computing Schur vectors if desired
632 * (Workspace: need N)
633 *
634  iwrk = itau
635  CALL shgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
636  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
637  $ work( iwrk ), lwork+1-iwrk, ierr )
638  IF( ierr.NE.0 ) THEN
639  IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
640  info = ierr
641  ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
642  info = ierr - n
643  ELSE
644  info = n + 1
645  END IF
646  GO TO 50
647  END IF
648 *
649 * Sort eigenvalues ALPHA/BETA and compute the reciprocal of
650 * condition number(s)
651 * (Workspace: If IJOB >= 1, need MAX( 8*(N+1), 2*SDIM*(N-SDIM) )
652 * otherwise, need 8*(N+1) )
653 *
654  IF( wantst ) THEN
655 *
656 * Undo scaling on eigenvalues before SELCTGing
657 *
658  IF( ilascl ) THEN
659  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n,
660  $ ierr )
661  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n,
662  $ ierr )
663  END IF
664  IF( ilbscl )
665  $ CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
666 *
667 * Select eigenvalues
668 *
669  DO 10 i = 1, n
670  bwork( i ) = selctg( alphar( i ), alphai( i ), beta( i ) )
671  10 CONTINUE
672 *
673 * Reorder eigenvalues, transform Generalized Schur vectors, and
674 * compute reciprocal condition numbers
675 *
676  CALL stgsen( ijob, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
677  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
678  $ sdim, pl, pr, dif, work( iwrk ), lwork-iwrk+1,
679  $ iwork, liwork, ierr )
680 *
681  IF( ijob.GE.1 )
682  $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
683  IF( ierr.EQ.-22 ) THEN
684 *
685 * not enough real workspace
686 *
687  info = -22
688  ELSE
689  IF( ijob.EQ.1 .OR. ijob.EQ.4 ) THEN
690  rconde( 1 ) = pl
691  rconde( 2 ) = pr
692  END IF
693  IF( ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
694  rcondv( 1 ) = dif( 1 )
695  rcondv( 2 ) = dif( 2 )
696  END IF
697  IF( ierr.EQ.1 )
698  $ info = n + 3
699  END IF
700 *
701  END IF
702 *
703 * Apply permutation to VSL and VSR
704 * (Workspace: none needed)
705 *
706  IF( ilvsl )
707  $ CALL sggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
708  $ work( iright ), n, vsl, ldvsl, ierr )
709 *
710  IF( ilvsr )
711  $ CALL sggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
712  $ work( iright ), n, vsr, ldvsr, ierr )
713 *
714 * Check if unscaling would cause over/underflow, if so, rescale
715 * (ALPHAR(I),ALPHAI(I),BETA(I)) so BETA(I) is on the order of
716 * B(I,I) and ALPHAR(I) and ALPHAI(I) are on the order of A(I,I)
717 *
718  IF( ilascl ) THEN
719  DO 20 i = 1, n
720  IF( alphai( i ).NE.zero ) THEN
721  IF( ( alphar( i ) / safmax ).GT.( anrmto / anrm ) .OR.
722  $ ( safmin / alphar( i ) ).GT.( anrm / anrmto ) )
723  $ THEN
724  work( 1 ) = abs( a( i, i ) / alphar( i ) )
725  beta( i ) = beta( i )*work( 1 )
726  alphar( i ) = alphar( i )*work( 1 )
727  alphai( i ) = alphai( i )*work( 1 )
728  ELSE IF( ( alphai( i ) / safmax ).GT.( anrmto / anrm )
729  $ .OR. ( safmin / alphai( i ) ).GT.( anrm / anrmto ) )
730  $ THEN
731  work( 1 ) = abs( a( i, i+1 ) / alphai( i ) )
732  beta( i ) = beta( i )*work( 1 )
733  alphar( i ) = alphar( i )*work( 1 )
734  alphai( i ) = alphai( i )*work( 1 )
735  END IF
736  END IF
737  20 CONTINUE
738  END IF
739 *
740  IF( ilbscl ) THEN
741  DO 25 i = 1, n
742  IF( alphai( i ).NE.zero ) THEN
743  IF( ( beta( i ) / safmax ).GT.( bnrmto / bnrm ) .OR.
744  $ ( safmin / beta( i ) ).GT.( bnrm / bnrmto ) ) THEN
745  work( 1 ) = abs( b( i, i ) / beta( i ) )
746  beta( i ) = beta( i )*work( 1 )
747  alphar( i ) = alphar( i )*work( 1 )
748  alphai( i ) = alphai( i )*work( 1 )
749  END IF
750  END IF
751  25 CONTINUE
752  END IF
753 *
754 * Undo scaling
755 *
756  IF( ilascl ) THEN
757  CALL slascl( 'H', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
758  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n, ierr )
759  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n, ierr )
760  END IF
761 *
762  IF( ilbscl ) THEN
763  CALL slascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
764  CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
765  END IF
766 *
767  IF( wantst ) THEN
768 *
769 * Check if reordering is correct
770 *
771  lastsl = .true.
772  lst2sl = .true.
773  sdim = 0
774  ip = 0
775  DO 40 i = 1, n
776  cursl = selctg( alphar( i ), alphai( i ), beta( i ) )
777  IF( alphai( i ).EQ.zero ) THEN
778  IF( cursl )
779  $ sdim = sdim + 1
780  ip = 0
781  IF( cursl .AND. .NOT.lastsl )
782  $ info = n + 2
783  ELSE
784  IF( ip.EQ.1 ) THEN
785 *
786 * Last eigenvalue of conjugate pair
787 *
788  cursl = cursl .OR. lastsl
789  lastsl = cursl
790  IF( cursl )
791  $ sdim = sdim + 2
792  ip = -1
793  IF( cursl .AND. .NOT.lst2sl )
794  $ info = n + 2
795  ELSE
796 *
797 * First eigenvalue of conjugate pair
798 *
799  ip = 1
800  END IF
801  END IF
802  lst2sl = lastsl
803  lastsl = cursl
804  40 CONTINUE
805 *
806  END IF
807 *
808  50 CONTINUE
809 *
810  work( 1 ) = maxwrk
811  iwork( 1 ) = liwmin
812 *
813  RETURN
814 *
815 * End of SGGESX
816 *
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:74
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:143
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:110
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:103
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine sggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
SGGBAK
Definition: sggbak.f:147
subroutine sggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
SGGBAL
Definition: sggbal.f:177
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:114
subroutine sgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
SGEQRF
Definition: sgeqrf.f:145
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:304
subroutine stgsen(IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, M, PL, PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO)
STGSEN
Definition: stgsen.f:451
subroutine sorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
SORGQR
Definition: sorgqr.f:128
subroutine sormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMQR
Definition: sormqr.f:168
subroutine sgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
SGGHRD
Definition: sgghrd.f:207
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:68
Here is the call graph for this function:
Here is the caller graph for this function: