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

◆ 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
!>  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 359 of file sggesx.f.

364*
365* -- LAPACK driver routine --
366* -- LAPACK is a software package provided by Univ. of Tennessee, --
367* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
368*
369* .. Scalar Arguments ..
370 CHARACTER JOBVSL, JOBVSR, SENSE, SORT
371 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LIWORK, LWORK, N,
372 $ SDIM
373* ..
374* .. Array Arguments ..
375 LOGICAL BWORK( * )
376 INTEGER IWORK( * )
377 REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
378 $ B( LDB, * ), BETA( * ), RCONDE( 2 ),
379 $ RCONDV( 2 ), VSL( LDVSL, * ), VSR( LDVSR, * ),
380 $ WORK( * )
381* ..
382* .. Function Arguments ..
383 LOGICAL SELCTG
384 EXTERNAL selctg
385* ..
386*
387* =====================================================================
388*
389* .. Parameters ..
390 REAL ZERO, ONE
391 parameter( zero = 0.0e+0, one = 1.0e+0 )
392* ..
393* .. Local Scalars ..
394 LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
395 $ LQUERY, LST2SL, WANTSB, WANTSE, WANTSN, WANTST,
396 $ WANTSV
397 INTEGER I, ICOLS, IERR, IHI, IJOB, IJOBVL, IJOBVR,
398 $ ILEFT, ILO, IP, IRIGHT, IROWS, ITAU, IWRK,
399 $ LIWMIN, LWRK, MAXWRK, MINWRK
400 REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PL,
401 $ PR, SAFMAX, SAFMIN, SMLNUM
402* ..
403* .. Local Arrays ..
404 REAL DIF( 2 )
405* ..
406* .. External Subroutines ..
407 EXTERNAL sgeqrf, sggbak, sggbal, sgghrd, shgeqz,
408 $ slacpy,
410* ..
411* .. External Functions ..
412 LOGICAL LSAME
413 INTEGER ILAENV
414 REAL SLAMCH, SLANGE, SROUNDUP_LWORK
415 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.
471 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
472 info = -3
473 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
474 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
475 info = -5
476 ELSE IF( n.LT.0 ) THEN
477 info = -6
478 ELSE IF( lda.LT.max( 1, n ) ) THEN
479 info = -8
480 ELSE IF( ldb.LT.max( 1, n ) ) THEN
481 info = -10
482 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
483 info = -16
484 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
485 info = -18
486 END IF
487*
488* Compute workspace
489* (Note: Comments in the code beginning "Workspace:" describe the
490* minimal amount of workspace needed at that point in the code,
491* as well as the preferred amount for good performance.
492* NB refers to the optimal block size for the immediately
493* following subroutine, as returned by ILAENV.)
494*
495 IF( info.EQ.0 ) THEN
496 IF( n.GT.0) THEN
497 minwrk = max( 8*n, 6*n + 16 )
498 maxwrk = minwrk - n +
499 $ n*ilaenv( 1, 'SGEQRF', ' ', n, 1, n, 0 )
500 maxwrk = max( maxwrk, minwrk - n +
501 $ n*ilaenv( 1, 'SORMQR', ' ', n, 1, n, -1 ) )
502 IF( ilvsl ) THEN
503 maxwrk = max( maxwrk, minwrk - n +
504 $ n*ilaenv( 1, 'SORGQR', ' ', n, 1, n, -1 ) )
505 END IF
506 lwrk = maxwrk
507 IF( ijob.GE.1 )
508 $ lwrk = max( lwrk, n*n/2 )
509 ELSE
510 minwrk = 1
511 maxwrk = 1
512 lwrk = 1
513 END IF
514 work( 1 ) = sroundup_lwork(lwrk)
515 IF( wantsn .OR. n.EQ.0 ) THEN
516 liwmin = 1
517 ELSE
518 liwmin = n + 6
519 END IF
520 iwork( 1 ) = liwmin
521*
522 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
523 info = -22
524 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
525 info = -24
526 END IF
527 END IF
528*
529 IF( info.NE.0 ) THEN
530 CALL xerbla( 'SGGESX', -info )
531 RETURN
532 ELSE IF (lquery) THEN
533 RETURN
534 END IF
535*
536* Quick return if possible
537*
538 IF( n.EQ.0 ) THEN
539 sdim = 0
540 RETURN
541 END IF
542*
543* Get machine constants
544*
545 eps = slamch( 'P' )
546 safmin = slamch( 'S' )
547 safmax = one / safmin
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,
666 $ ierr )
667*
668* Select eigenvalues
669*
670 DO 10 i = 1, n
671 bwork( i ) = selctg( alphar( i ), alphai( i ),
672 $ beta( i ) )
673 10 CONTINUE
674*
675* Reorder eigenvalues, transform Generalized Schur vectors, and
676* compute reciprocal condition numbers
677*
678 CALL stgsen( ijob, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
679 $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
680 $ sdim, pl, pr, dif, work( iwrk ), lwork-iwrk+1,
681 $ iwork, liwork, ierr )
682*
683 IF( ijob.GE.1 )
684 $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
685 IF( ierr.EQ.-22 ) THEN
686*
687* not enough real workspace
688*
689 info = -22
690 ELSE
691 IF( ijob.EQ.1 .OR. ijob.EQ.4 ) THEN
692 rconde( 1 ) = pl
693 rconde( 2 ) = pr
694 END IF
695 IF( ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
696 rcondv( 1 ) = dif( 1 )
697 rcondv( 2 ) = dif( 2 )
698 END IF
699 IF( ierr.EQ.1 )
700 $ info = n + 3
701 END IF
702*
703 END IF
704*
705* Apply permutation to VSL and VSR
706* (Workspace: none needed)
707*
708 IF( ilvsl )
709 $ CALL sggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
710 $ work( iright ), n, vsl, ldvsl, ierr )
711*
712 IF( ilvsr )
713 $ CALL sggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
714 $ work( iright ), n, vsr, ldvsr, ierr )
715*
716* Check if unscaling would cause over/underflow, if so, rescale
717* (ALPHAR(I),ALPHAI(I),BETA(I)) so BETA(I) is on the order of
718* B(I,I) and ALPHAR(I) and ALPHAI(I) are on the order of A(I,I)
719*
720 IF( ilascl ) THEN
721 DO 20 i = 1, n
722 IF( alphai( i ).NE.zero ) THEN
723 IF( ( alphar( i ) / safmax ).GT.( anrmto / anrm ) .OR.
724 $ ( safmin / alphar( i ) ).GT.( anrm / anrmto ) )
725 $ THEN
726 work( 1 ) = abs( a( i, i ) / alphar( i ) )
727 beta( i ) = beta( i )*work( 1 )
728 alphar( i ) = alphar( i )*work( 1 )
729 alphai( i ) = alphai( i )*work( 1 )
730 ELSE IF( ( alphai( i ) / safmax ).GT.( anrmto / anrm )
731 $ .OR. ( safmin / alphai( i ) ).GT.( anrm / anrmto ) )
732 $ THEN
733 work( 1 ) = abs( a( i, i+1 ) / alphai( i ) )
734 beta( i ) = beta( i )*work( 1 )
735 alphar( i ) = alphar( i )*work( 1 )
736 alphai( i ) = alphai( i )*work( 1 )
737 END IF
738 END IF
739 20 CONTINUE
740 END IF
741*
742 IF( ilbscl ) THEN
743 DO 25 i = 1, n
744 IF( alphai( i ).NE.zero ) THEN
745 IF( ( beta( i ) / safmax ).GT.( bnrmto / bnrm ) .OR.
746 $ ( safmin / beta( i ) ).GT.( bnrm / bnrmto ) ) THEN
747 work( 1 ) = abs( b( i, i ) / beta( i ) )
748 beta( i ) = beta( i )*work( 1 )
749 alphar( i ) = alphar( i )*work( 1 )
750 alphai( i ) = alphai( i )*work( 1 )
751 END IF
752 END IF
753 25 CONTINUE
754 END IF
755*
756* Undo scaling
757*
758 IF( ilascl ) THEN
759 CALL slascl( 'H', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
760 CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n,
761 $ ierr )
762 CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n,
763 $ ierr )
764 END IF
765*
766 IF( ilbscl ) THEN
767 CALL slascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
768 CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
769 END IF
770*
771 IF( wantst ) THEN
772*
773* Check if reordering is correct
774*
775 lastsl = .true.
776 lst2sl = .true.
777 sdim = 0
778 ip = 0
779 DO 40 i = 1, n
780 cursl = selctg( alphar( i ), alphai( i ), beta( i ) )
781 IF( alphai( i ).EQ.zero ) THEN
782 IF( cursl )
783 $ sdim = sdim + 1
784 ip = 0
785 IF( cursl .AND. .NOT.lastsl )
786 $ info = n + 2
787 ELSE
788 IF( ip.EQ.1 ) THEN
789*
790* Last eigenvalue of conjugate pair
791*
792 cursl = cursl .OR. lastsl
793 lastsl = cursl
794 IF( cursl )
795 $ sdim = sdim + 2
796 ip = -1
797 IF( cursl .AND. .NOT.lst2sl )
798 $ info = n + 2
799 ELSE
800*
801* First eigenvalue of conjugate pair
802*
803 ip = 1
804 END IF
805 END IF
806 lst2sl = lastsl
807 lastsl = cursl
808 40 CONTINUE
809*
810 END IF
811*
812 50 CONTINUE
813*
814 work( 1 ) = sroundup_lwork(maxwrk)
815 iwork( 1 ) = liwmin
816*
817 RETURN
818*
819* End of SGGESX
820*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgeqrf(m, n, a, lda, tau, work, lwork, info)
SGEQRF
Definition sgeqrf.f:144
subroutine sggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
SGGBAK
Definition sggbak.f:146
subroutine sggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
SGGBAL
Definition sggbal.f:175
subroutine sgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
SGGHRD
Definition sgghrd.f:206
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:303
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
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
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:108
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
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:450
subroutine sorgqr(m, n, k, a, lda, tau, work, lwork, info)
SORGQR
Definition sorgqr.f:126
subroutine sormqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
SORMQR
Definition sormqr.f:166
Here is the call graph for this function:
Here is the caller graph for this function: