387 SUBROUTINE dggevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
388 $ ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO,
389 $ IHI, LSCALE, RSCALE, ABNRM, BBNRM, RCONDE,
390 $ RCONDV, WORK, LWORK, IWORK, BWORK, INFO )
397 CHARACTER BALANC, JOBVL, JOBVR, SENSE
398 INTEGER IHI, ILO, INFO, LDA, LDB, LDVL, LDVR, LWORK, N
399 DOUBLE PRECISION ABNRM, BBNRM
404 DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
405 $ b( ldb, * ), beta( * ), lscale( * ),
406 $ rconde( * ), rcondv( * ), rscale( * ),
407 $ vl( ldvl, * ), vr( ldvr, * ), work( * )
413 DOUBLE PRECISION ZERO, ONE
414 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
417 LOGICAL ILASCL, ILBSCL, ILV, ILVL, ILVR, LQUERY, NOSCL,
418 $ PAIR, WANTSB, WANTSE, WANTSN, WANTSV
420 INTEGER I, ICOLS, IERR, IJOBVL, IJOBVR, IN, IROWS,
421 $ ITAU, IWRK, IWRK1, J, JC, JR, M, MAXWRK,
423 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
437 DOUBLE PRECISION DLAMCH, DLANGE
438 EXTERNAL lsame, ilaenv, dlamch, dlange
441 INTRINSIC abs, max, sqrt
447 IF( lsame( jobvl,
'N' ) )
THEN
450 ELSE IF( lsame( jobvl,
'V' ) )
THEN
458 IF( lsame( jobvr,
'N' ) )
THEN
461 ELSE IF( lsame( jobvr,
'V' ) )
THEN
470 noscl = lsame( balanc,
'N' ) .OR. lsame( balanc,
'P' )
471 wantsn = lsame( sense,
'N' )
472 wantse = lsame( sense,
'E' )
473 wantsv = lsame( sense,
'V' )
474 wantsb = lsame( sense,
'B' )
479 lquery = ( lwork.EQ.-1 )
480 IF( .NOT.( lsame( balanc,
'N' ) .OR. lsame( balanc,
481 $
'S' ) .OR. lsame( balanc,
'P' ) .OR. lsame( balanc,
'B' ) ) )
484 ELSE IF( ijobvl.LE.0 )
THEN
486 ELSE IF( ijobvr.LE.0 )
THEN
488 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsb .OR. wantsv ) )
491 ELSE IF( n.LT.0 )
THEN
493 ELSE IF( lda.LT.max( 1, n ) )
THEN
495 ELSE IF( ldb.LT.max( 1, n ) )
THEN
497 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) )
THEN
499 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) )
THEN
516 IF( noscl .AND. .NOT.ilv )
THEN
521 IF( wantse .OR. wantsb )
THEN
524 IF( wantsv .OR. wantsb )
THEN
525 minwrk = max( minwrk, 2*n*( n + 4 ) + 16 )
528 maxwrk = max( maxwrk,
529 $ n + n*ilaenv( 1,
'DGEQRF',
' ', n, 1, n, 0 ) )
530 maxwrk = max( maxwrk,
531 $ n + n*ilaenv( 1,
'DORMQR',
' ', n, 1, n, 0 ) )
533 maxwrk = max( maxwrk, n +
534 $ n*ilaenv( 1,
'DORGQR',
' ', n, 1, n, 0 ) )
539 IF( lwork.LT.minwrk .AND. .NOT.lquery )
THEN
545 CALL xerbla(
'DGGEVX', -info )
547 ELSE IF( lquery )
THEN
560 smlnum = dlamch(
'S' )
561 bignum = one / smlnum
562 CALL dlabad( smlnum, bignum )
563 smlnum = sqrt( smlnum ) / eps
564 bignum = one / smlnum
568 anrm = dlange(
'M', n, n, a, lda, work )
570 IF( anrm.GT.zero .AND. anrm.LT.smlnum )
THEN
573 ELSE IF( anrm.GT.bignum )
THEN
578 $
CALL dlascl(
'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
582 bnrm = dlange(
'M', n, n, b, ldb, work )
584 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum )
THEN
587 ELSE IF( bnrm.GT.bignum )
THEN
592 $
CALL dlascl(
'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
597 CALL dggbal( balanc, n, a, lda, b, ldb, ilo, ihi, lscale, rscale,
602 abnrm = dlange(
'1', n, n, a, lda, work( 1 ) )
605 CALL dlascl(
'G', 0, 0, anrmto, anrm, 1, 1, work( 1 ), 1,
610 bbnrm = dlange(
'1', n, n, b, ldb, work( 1 ) )
613 CALL dlascl(
'G', 0, 0, bnrmto, bnrm, 1, 1, work( 1 ), 1,
621 irows = ihi + 1 - ilo
622 IF( ilv .OR. .NOT.wantsn )
THEN
629 CALL dgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
630 $ work( iwrk ), lwork+1-iwrk, ierr )
635 CALL dormqr(
'L',
'T', irows, icols, irows, b( ilo, ilo ), ldb,
636 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
637 $ lwork+1-iwrk, ierr )
643 CALL dlaset(
'Full', n, n, zero, one, vl, ldvl )
644 IF( irows.GT.1 )
THEN
645 CALL dlacpy(
'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
646 $ vl( ilo+1, ilo ), ldvl )
648 CALL dorgqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
649 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
653 $
CALL dlaset(
'Full', n, n, zero, one, vr, ldvr )
658 IF( ilv .OR. .NOT.wantsn )
THEN
662 CALL dgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
663 $ ldvl, vr, ldvr, ierr )
665 CALL dgghrd(
'N',
'N', irows, 1, irows, a( ilo, ilo ), lda,
666 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
673 IF( ilv .OR. .NOT.wantsn )
THEN
679 CALL dhgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
680 $ alphar, alphai, beta, vl, ldvl, vr, ldvr, work,
683 IF( ierr.GT.0 .AND. ierr.LE.n )
THEN
685 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n )
THEN
698 IF( ilv .OR. .NOT.wantsn )
THEN
710 CALL dtgevc( chtemp,
'B', ldumma, n, a, lda, b, ldb, vl,
711 $ ldvl, vr, ldvr, n, in, work, ierr )
718 IF( .NOT.wantsn )
THEN
737 IF( a( i+1, i ).NE.zero )
THEN
748 ELSE IF( mm.EQ.2 )
THEN
750 bwork( i+1 ) = .true.
759 IF( wantse .OR. wantsb )
THEN
760 CALL dtgevc(
'B',
'S', bwork, n, a, lda, b, ldb,
761 $ work( 1 ), n, work( iwrk ), n, mm, m,
762 $ work( iwrk1 ), ierr )
769 CALL dtgsna( sense,
'S', bwork, n, a, lda, b, ldb,
770 $ work( 1 ), n, work( iwrk ), n, rconde( i ),
771 $ rcondv( i ), mm, m, work( iwrk1 ),
772 $ lwork-iwrk1+1, iwork, ierr )
782 CALL dggbak( balanc,
'L', n, ilo, ihi, lscale, rscale, n, vl,
786 IF( alphai( jc ).LT.zero )
789 IF( alphai( jc ).EQ.zero )
THEN
791 temp = max( temp, abs( vl( jr, jc ) ) )
795 temp = max( temp, abs( vl( jr, jc ) )+
796 $ abs( vl( jr, jc+1 ) ) )
802 IF( alphai( jc ).EQ.zero )
THEN
804 vl( jr, jc ) = vl( jr, jc )*temp
808 vl( jr, jc ) = vl( jr, jc )*temp
809 vl( jr, jc+1 ) = vl( jr, jc+1 )*temp
815 CALL dggbak( balanc,
'R', n, ilo, ihi, lscale, rscale, n, vr,
818 IF( alphai( jc ).LT.zero )
821 IF( alphai( jc ).EQ.zero )
THEN
823 temp = max( temp, abs( vr( jr, jc ) ) )
827 temp = max( temp, abs( vr( jr, jc ) )+
828 $ abs( vr( jr, jc+1 ) ) )
834 IF( alphai( jc ).EQ.zero )
THEN
836 vr( jr, jc ) = vr( jr, jc )*temp
840 vr( jr, jc ) = vr( jr, jc )*temp
841 vr( jr, jc+1 ) = vr( jr, jc+1 )*temp
852 CALL dlascl(
'G', 0, 0, anrmto, anrm, n, 1, alphar, n, ierr )
853 CALL dlascl(
'G', 0, 0, anrmto, anrm, n, 1, alphai, n, ierr )
857 CALL dlascl(
'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
subroutine dlabad(SMALL, LARGE)
DLABAD
subroutine dlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine dggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
DGGBAK
subroutine dggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
DGGBAL
subroutine dhgeqz(JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, INFO)
DHGEQZ
subroutine dgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
DGEQRF
subroutine dtgevc(SIDE, HOWMNY, SELECT, N, S, LDS, P, LDP, VL, LDVL, VR, LDVR, MM, M, WORK, INFO)
DTGEVC
subroutine dggevx(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, BWORK, INFO)
DGGEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
subroutine dorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGQR
subroutine dgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
DGGHRD
subroutine dormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMQR
subroutine dtgsna(JOB, HOWMNY, SELECT, N, A, LDA, B, LDB, VL, LDVL, VR, LDVR, S, DIF, MM, M, WORK, LWORK, IWORK, INFO)
DTGSNA