304 SUBROUTINE dgeevx( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, WI,
305 $ vl, ldvl, vr, ldvr, ilo, ihi, scale, abnrm,
306 $ rconde, rcondv, work, lwork, iwork, info )
315 CHARACTER BALANC, JOBVL, JOBVR, SENSE
316 INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
317 DOUBLE PRECISION ABNRM
321 DOUBLE PRECISION A( lda, * ), RCONDE( * ), RCONDV( * ),
322 $ scale( * ), vl( ldvl, * ), vr( ldvr, * ),
323 $ wi( * ), work( * ), wr( * )
329 DOUBLE PRECISION ZERO, ONE
330 parameter ( zero = 0.0d0, one = 1.0d0 )
333 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
336 INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
337 $ lwork_trevc, maxwrk, minwrk, nout
338 DOUBLE PRECISION ANRM, BIGNUM, CS, CSCALE, EPS, R, SCL, SMLNUM,
343 DOUBLE PRECISION DUM( 1 )
352 INTEGER IDAMAX, ILAENV
353 DOUBLE PRECISION DLAMCH, DLANGE, DLAPY2, DNRM2
354 EXTERNAL lsame, idamax, ilaenv, dlamch, dlange, dlapy2,
365 lquery = ( lwork.EQ.-1 )
366 wantvl = lsame( jobvl,
'V' )
367 wantvr = lsame( jobvr,
'V' )
368 wntsnn = lsame( sense,
'N' )
369 wntsne = lsame( sense,
'E' )
370 wntsnv = lsame( sense,
'V' )
371 wntsnb = lsame( sense,
'B' )
372 IF( .NOT.( lsame( balanc,
'N' ) .OR. lsame( balanc,
'S' )
373 $ .OR. lsame( balanc,
'P' ) .OR. lsame( balanc,
'B' ) ) )
376 ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl,
'N' ) ) )
THEN
378 ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr,
'N' ) ) )
THEN
380 ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
381 $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
384 ELSE IF( n.LT.0 )
THEN
386 ELSE IF( lda.LT.max( 1, n ) )
THEN
388 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) )
THEN
390 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) )
THEN
409 maxwrk = n + n*ilaenv( 1,
'DGEHRD',
' ', n, 1, n, 0 )
412 CALL dtrevc3(
'L',
'B',
SELECT, n, a, lda,
413 $ vl, ldvl, vr, ldvr,
414 $ n, nout, work, -1, ierr )
415 lwork_trevc = int( work(1) )
416 maxwrk = max( maxwrk, n + lwork_trevc )
417 CALL dhseqr(
'S',
'V', n, 1, n, a, lda, wr, wi, vl, ldvl,
419 ELSE IF( wantvr )
THEN
420 CALL dtrevc3(
'R',
'B',
SELECT, n, a, lda,
421 $ vl, ldvl, vr, ldvr,
422 $ n, nout, work, -1, ierr )
423 lwork_trevc = int( work(1) )
424 maxwrk = max( maxwrk, n + lwork_trevc )
425 CALL dhseqr(
'S',
'V', n, 1, n, a, lda, wr, wi, vr, ldvr,
429 CALL dhseqr(
'E',
'N', n, 1, n, a, lda, wr, wi, vr,
430 $ ldvr, work, -1, info )
432 CALL dhseqr(
'S',
'N', n, 1, n, a, lda, wr, wi, vr,
433 $ ldvr, work, -1, info )
436 hswork = int( work(1) )
438 IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) )
THEN
441 $ minwrk = max( minwrk, n*n+6*n )
442 maxwrk = max( maxwrk, hswork )
444 $ maxwrk = max( maxwrk, n*n + 6*n )
447 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
448 $ minwrk = max( minwrk, n*n + 6*n )
449 maxwrk = max( maxwrk, hswork )
450 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
'DORGHR',
451 $
' ', n, 1, n, -1 ) )
452 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
453 $ maxwrk = max( maxwrk, n*n + 6*n )
454 maxwrk = max( maxwrk, 3*n )
456 maxwrk = max( maxwrk, minwrk )
460 IF( lwork.LT.minwrk .AND. .NOT.lquery )
THEN
466 CALL xerbla(
'DGEEVX', -info )
468 ELSE IF( lquery )
THEN
480 smlnum = dlamch(
'S' )
481 bignum = one / smlnum
482 CALL dlabad( smlnum, bignum )
483 smlnum = sqrt( smlnum ) / eps
484 bignum = one / smlnum
489 anrm = dlange(
'M', n, n, a, lda, dum )
491 IF( anrm.GT.zero .AND. anrm.LT.smlnum )
THEN
494 ELSE IF( anrm.GT.bignum )
THEN
499 $
CALL dlascl(
'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
503 CALL dgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
504 abnrm = dlange(
'1', n, n, a, lda, dum )
507 CALL dlascl(
'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
516 CALL dgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
517 $ lwork-iwrk+1, ierr )
525 CALL dlacpy(
'L', n, n, a, lda, vl, ldvl )
530 CALL dorghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
531 $ lwork-iwrk+1, ierr )
537 CALL dhseqr(
'S',
'V', n, ilo, ihi, a, lda, wr, wi, vl, ldvl,
538 $ work( iwrk ), lwork-iwrk+1, info )
546 CALL dlacpy(
'F', n, n, vl, ldvl, vr, ldvr )
549 ELSE IF( wantvr )
THEN
555 CALL dlacpy(
'L', n, n, a, lda, vr, ldvr )
560 CALL dorghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
561 $ lwork-iwrk+1, ierr )
567 CALL dhseqr(
'S',
'V', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
568 $ work( iwrk ), lwork-iwrk+1, info )
584 CALL dhseqr( job,
'N', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
585 $ work( iwrk ), lwork-iwrk+1, info )
593 IF( wantvl .OR. wantvr )
THEN
598 CALL dtrevc3( side,
'B',
SELECT, n, a, lda, vl, ldvl, vr, ldvr,
599 $ n, nout, work( iwrk ), lwork-iwrk+1, ierr )
605 IF( .NOT.wntsnn )
THEN
606 CALL dtrsna( sense,
'A',
SELECT, n, a, lda, vl, ldvl, vr, ldvr,
607 $ rconde, rcondv, n, nout, work( iwrk ), n, iwork,
615 CALL dgebak( balanc,
'L', n, ilo, ihi, scale, n, vl, ldvl,
621 IF( wi( i ).EQ.zero )
THEN
622 scl = one / dnrm2( n, vl( 1, i ), 1 )
623 CALL dscal( n, scl, vl( 1, i ), 1 )
624 ELSE IF( wi( i ).GT.zero )
THEN
625 scl = one / dlapy2( dnrm2( n, vl( 1, i ), 1 ),
626 $ dnrm2( n, vl( 1, i+1 ), 1 ) )
627 CALL dscal( n, scl, vl( 1, i ), 1 )
628 CALL dscal( n, scl, vl( 1, i+1 ), 1 )
630 work( k ) = vl( k, i )**2 + vl( k, i+1 )**2
632 k = idamax( n, work, 1 )
633 CALL dlartg( vl( k, i ), vl( k, i+1 ), cs, sn, r )
634 CALL drot( n, vl( 1, i ), 1, vl( 1, i+1 ), 1, cs, sn )
644 CALL dgebak( balanc,
'R', n, ilo, ihi, scale, n, vr, ldvr,
650 IF( wi( i ).EQ.zero )
THEN
651 scl = one / dnrm2( n, vr( 1, i ), 1 )
652 CALL dscal( n, scl, vr( 1, i ), 1 )
653 ELSE IF( wi( i ).GT.zero )
THEN
654 scl = one / dlapy2( dnrm2( n, vr( 1, i ), 1 ),
655 $ dnrm2( n, vr( 1, i+1 ), 1 ) )
656 CALL dscal( n, scl, vr( 1, i ), 1 )
657 CALL dscal( n, scl, vr( 1, i+1 ), 1 )
659 work( k ) = vr( k, i )**2 + vr( k, i+1 )**2
661 k = idamax( n, work, 1 )
662 CALL dlartg( vr( k, i ), vr( k, i+1 ), cs, sn, r )
663 CALL drot( n, vr( 1, i ), 1, vr( 1, i+1 ), 1, cs, sn )
673 CALL dlascl(
'G', 0, 0, cscale, anrm, n-info, 1, wr( info+1 ),
674 $ max( n-info, 1 ), ierr )
675 CALL dlascl(
'G', 0, 0, cscale, anrm, n-info, 1, wi( info+1 ),
676 $ max( n-info, 1 ), ierr )
678 IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
679 $
CALL dlascl(
'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
682 CALL dlascl(
'G', 0, 0, cscale, anrm, ilo-1, 1, wr, n,
684 CALL dlascl(
'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
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 dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
subroutine dtrevc3(SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, LWORK, INFO)
DTREVC3
subroutine drot(N, DX, INCX, DY, INCY, C, S)
DROT
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
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 dgebak(JOB, SIDE, N, ILO, IHI, SCALE, M, V, LDV, INFO)
DGEBAK
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine dlabad(SMALL, LARGE)
DLABAD
subroutine dgebal(JOB, N, A, LDA, ILO, IHI, SCALE, INFO)
DGEBAL
subroutine dtrsna(JOB, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, S, SEP, MM, M, WORK, LDWORK, IWORK, INFO)
DTRSNA
subroutine dorghr(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DORGHR
subroutine dscal(N, DA, DX, INCX)
DSCAL
subroutine dhseqr(JOB, COMPZ, N, ILO, IHI, H, LDH, WR, WI, Z, LDZ, WORK, LWORK, INFO)
DHSEQR
subroutine dlartg(F, G, CS, SN, R)
DLARTG generates a plane rotation with real cosine and real sine.