303 SUBROUTINE dgeevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, WI,
304 $ VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM,
305 $ RCONDE, RCONDV, WORK, LWORK, IWORK, INFO )
313 CHARACTER BALANC, JOBVL, JOBVR, SENSE
314 INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
315 DOUBLE PRECISION ABNRM
319 DOUBLE PRECISION A( LDA, * ), RCONDE( * ), RCONDV( * ),
320 $ scale( * ), vl( ldvl, * ), vr( ldvr, * ),
321 $ wi( * ), work( * ), wr( * )
327 DOUBLE PRECISION ZERO, ONE
328 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
331 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
334 INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
335 $ lwork_trevc, maxwrk, minwrk, nout
336 DOUBLE PRECISION ANRM, BIGNUM, CS, CSCALE, EPS, R, SCL, SMLNUM,
341 DOUBLE PRECISION DUM( 1 )
350 INTEGER IDAMAX, ILAENV
351 DOUBLE PRECISION DLAMCH, DLANGE, DLAPY2, DNRM2
352 EXTERNAL lsame, idamax, ilaenv, dlamch, dlange, dlapy2,
363 lquery = ( lwork.EQ.-1 )
364 wantvl = lsame( jobvl,
'V' )
365 wantvr = lsame( jobvr,
'V' )
366 wntsnn = lsame( sense,
'N' )
367 wntsne = lsame( sense,
'E' )
368 wntsnv = lsame( sense,
'V' )
369 wntsnb = lsame( sense,
'B' )
370 IF( .NOT.( lsame( balanc,
'N' ) .OR. lsame( balanc,
'S' )
371 $ .OR. lsame( balanc,
'P' ) .OR. lsame( balanc,
'B' ) ) )
374 ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl,
'N' ) ) )
THEN
376 ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr,
'N' ) ) )
THEN
378 ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
379 $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
382 ELSE IF( n.LT.0 )
THEN
384 ELSE IF( lda.LT.max( 1, n ) )
THEN
386 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) )
THEN
388 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) )
THEN
407 maxwrk = n + n*ilaenv( 1,
'DGEHRD',
' ', n, 1, n, 0 )
410 CALL dtrevc3(
'L',
'B',
SELECT, n, a, lda,
411 $ vl, ldvl, vr, ldvr,
412 $ n, nout, work, -1, ierr )
413 lwork_trevc = int( work(1) )
414 maxwrk = max( maxwrk, n + lwork_trevc )
415 CALL dhseqr(
'S',
'V', n, 1, n, a, lda, wr, wi, vl, ldvl,
417 ELSE IF( wantvr )
THEN
418 CALL dtrevc3(
'R',
'B',
SELECT, n, a, lda,
419 $ vl, ldvl, vr, ldvr,
420 $ n, nout, work, -1, ierr )
421 lwork_trevc = int( work(1) )
422 maxwrk = max( maxwrk, n + lwork_trevc )
423 CALL dhseqr(
'S',
'V', n, 1, n, a, lda, wr, wi, vr, ldvr,
427 CALL dhseqr(
'E',
'N', n, 1, n, a, lda, wr, wi, vr,
428 $ ldvr, work, -1, info )
430 CALL dhseqr(
'S',
'N', n, 1, n, a, lda, wr, wi, vr,
431 $ ldvr, work, -1, info )
434 hswork = int( work(1) )
436 IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) )
THEN
439 $ minwrk = max( minwrk, n*n+6*n )
440 maxwrk = max( maxwrk, hswork )
442 $ maxwrk = max( maxwrk, n*n + 6*n )
445 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
446 $ minwrk = max( minwrk, n*n + 6*n )
447 maxwrk = max( maxwrk, hswork )
448 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
'DORGHR',
449 $
' ', n, 1, n, -1 ) )
450 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
451 $ maxwrk = max( maxwrk, n*n + 6*n )
452 maxwrk = max( maxwrk, 3*n )
454 maxwrk = max( maxwrk, minwrk )
458 IF( lwork.LT.minwrk .AND. .NOT.lquery )
THEN
464 CALL xerbla(
'DGEEVX', -info )
466 ELSE IF( lquery )
THEN
478 smlnum = dlamch(
'S' )
479 bignum = one / smlnum
480 CALL dlabad( smlnum, bignum )
481 smlnum = sqrt( smlnum ) / eps
482 bignum = one / smlnum
487 anrm = dlange(
'M', n, n, a, lda, dum )
489 IF( anrm.GT.zero .AND. anrm.LT.smlnum )
THEN
492 ELSE IF( anrm.GT.bignum )
THEN
497 $
CALL dlascl(
'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
501 CALL dgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
502 abnrm = dlange(
'1', n, n, a, lda, dum )
505 CALL dlascl(
'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
514 CALL dgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
515 $ lwork-iwrk+1, ierr )
523 CALL dlacpy(
'L', n, n, a, lda, vl, ldvl )
528 CALL dorghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
529 $ lwork-iwrk+1, ierr )
535 CALL dhseqr(
'S',
'V', n, ilo, ihi, a, lda, wr, wi, vl, ldvl,
536 $ work( iwrk ), lwork-iwrk+1, info )
544 CALL dlacpy(
'F', n, n, vl, ldvl, vr, ldvr )
547 ELSE IF( wantvr )
THEN
553 CALL dlacpy(
'L', n, n, a, lda, vr, ldvr )
558 CALL dorghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
559 $ lwork-iwrk+1, ierr )
565 CALL dhseqr(
'S',
'V', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
566 $ work( iwrk ), lwork-iwrk+1, info )
582 CALL dhseqr( job,
'N', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
583 $ work( iwrk ), lwork-iwrk+1, info )
591 IF( wantvl .OR. wantvr )
THEN
596 CALL dtrevc3( side,
'B',
SELECT, n, a, lda, vl, ldvl, vr, ldvr,
597 $ n, nout, work( iwrk ), lwork-iwrk+1, ierr )
603 IF( .NOT.wntsnn )
THEN
604 CALL dtrsna( sense,
'A',
SELECT, n, a, lda, vl, ldvl, vr, ldvr,
605 $ rconde, rcondv, n, nout, work( iwrk ), n, iwork,
613 CALL dgebak( balanc,
'L', n, ilo, ihi, scale, n, vl, ldvl,
619 IF( wi( i ).EQ.zero )
THEN
620 scl = one / dnrm2( n, vl( 1, i ), 1 )
621 CALL dscal( n, scl, vl( 1, i ), 1 )
622 ELSE IF( wi( i ).GT.zero )
THEN
623 scl = one / dlapy2( dnrm2( n, vl( 1, i ), 1 ),
624 $ dnrm2( n, vl( 1, i+1 ), 1 ) )
625 CALL dscal( n, scl, vl( 1, i ), 1 )
626 CALL dscal( n, scl, vl( 1, i+1 ), 1 )
628 work( k ) = vl( k, i )**2 + vl( k, i+1 )**2
630 k = idamax( n, work, 1 )
631 CALL dlartg( vl( k, i ), vl( k, i+1 ), cs, sn, r )
632 CALL drot( n, vl( 1, i ), 1, vl( 1, i+1 ), 1, cs, sn )
642 CALL dgebak( balanc,
'R', n, ilo, ihi, scale, n, vr, ldvr,
648 IF( wi( i ).EQ.zero )
THEN
649 scl = one / dnrm2( n, vr( 1, i ), 1 )
650 CALL dscal( n, scl, vr( 1, i ), 1 )
651 ELSE IF( wi( i ).GT.zero )
THEN
652 scl = one / dlapy2( dnrm2( n, vr( 1, i ), 1 ),
653 $ dnrm2( n, vr( 1, i+1 ), 1 ) )
654 CALL dscal( n, scl, vr( 1, i ), 1 )
655 CALL dscal( n, scl, vr( 1, i+1 ), 1 )
657 work( k ) = vr( k, i )**2 + vr( k, i+1 )**2
659 k = idamax( n, work, 1 )
660 CALL dlartg( vr( k, i ), vr( k, i+1 ), cs, sn, r )
661 CALL drot( n, vr( 1, i ), 1, vr( 1, i+1 ), 1, cs, sn )
671 CALL dlascl(
'G', 0, 0, cscale, anrm, n-info, 1, wr( info+1 ),
672 $ max( n-info, 1 ), ierr )
673 CALL dlascl(
'G', 0, 0, cscale, anrm, n-info, 1, wi( info+1 ),
674 $ max( n-info, 1 ), ierr )
676 IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
677 $
CALL dlascl(
'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
680 CALL dlascl(
'G', 0, 0, cscale, anrm, ilo-1, 1, wr, n,
682 CALL dlascl(
'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
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 dlartg(f, g, c, s, r)
DLARTG generates a plane rotation with real cosine and real sine.
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine drot(N, DX, INCX, DY, INCY, C, S)
DROT
subroutine dscal(N, DA, DX, INCX)
DSCAL
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
subroutine dgebal(JOB, N, A, LDA, ILO, IHI, SCALE, INFO)
DGEBAL
subroutine dgebak(JOB, SIDE, N, ILO, IHI, SCALE, M, V, LDV, INFO)
DGEBAK
subroutine dgeevx(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, WI, VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, INFO)
DGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
subroutine dtrsna(JOB, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, S, SEP, MM, M, WORK, LDWORK, IWORK, INFO)
DTRSNA
subroutine dhseqr(JOB, COMPZ, N, ILO, IHI, H, LDH, WR, WI, Z, LDZ, WORK, LWORK, INFO)
DHSEQR
subroutine dorghr(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DORGHR
subroutine dtrevc3(SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, LWORK, INFO)
DTREVC3