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

◆ zggevx()

subroutine zggevx ( character  balanc,
character  jobvl,
character  jobvr,
character  sense,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( ldb, * )  b,
integer  ldb,
complex*16, dimension( * )  alpha,
complex*16, dimension( * )  beta,
complex*16, dimension( ldvl, * )  vl,
integer  ldvl,
complex*16, dimension( ldvr, * )  vr,
integer  ldvr,
integer  ilo,
integer  ihi,
double precision, dimension( * )  lscale,
double precision, dimension( * )  rscale,
double precision  abnrm,
double precision  bbnrm,
double precision, dimension( * )  rconde,
double precision, dimension( * )  rcondv,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
logical, dimension( * )  bwork,
integer  info 
)

ZGGEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices

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

Purpose:
 ZGGEVX computes for a pair of N-by-N complex nonsymmetric matrices
 (A,B) the generalized eigenvalues, and optionally, the left and/or
 right generalized eigenvectors.

 Optionally, it also computes a balancing transformation to improve
 the conditioning of the eigenvalues and eigenvectors (ILO, IHI,
 LSCALE, RSCALE, ABNRM, and BBNRM), reciprocal condition numbers for
 the eigenvalues (RCONDE), and reciprocal condition numbers for the
 right eigenvectors (RCONDV).

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

 The right eigenvector v(j) corresponding to the eigenvalue lambda(j)
 of (A,B) satisfies
                  A * v(j) = lambda(j) * B * v(j) .
 The left eigenvector u(j) corresponding to the eigenvalue lambda(j)
 of (A,B) satisfies
                  u(j)**H * A  = lambda(j) * u(j)**H * B.
 where u(j)**H is the conjugate-transpose of u(j).
Parameters
[in]BALANC
          BALANC is CHARACTER*1
          Specifies the balance option to be performed:
          = 'N':  do not diagonally scale or permute;
          = 'P':  permute only;
          = 'S':  scale only;
          = 'B':  both permute and scale.
          Computed reciprocal condition numbers will be for the
          matrices after permuting and/or balancing. Permuting does
          not change condition numbers (in exact arithmetic), but
          balancing does.
[in]JOBVL
          JOBVL is CHARACTER*1
          = 'N':  do not compute the left generalized eigenvectors;
          = 'V':  compute the left generalized eigenvectors.
[in]JOBVR
          JOBVR is CHARACTER*1
          = 'N':  do not compute the right generalized eigenvectors;
          = 'V':  compute the right generalized eigenvectors.
[in]SENSE
          SENSE is CHARACTER*1
          Determines which reciprocal condition numbers are computed.
          = 'N': none are computed;
          = 'E': computed for eigenvalues only;
          = 'V': computed for eigenvectors only;
          = 'B': computed for eigenvalues and eigenvectors.
[in]N
          N is INTEGER
          The order of the matrices A, B, VL, and VR.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA, N)
          On entry, the matrix A in the pair (A,B).
          On exit, A has been overwritten. If JOBVL='V' or JOBVR='V'
          or both, then A contains the first part of the complex Schur
          form of the "balanced" versions of the input A and B.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB, N)
          On entry, the matrix B in the pair (A,B).
          On exit, B has been overwritten. If JOBVL='V' or JOBVR='V'
          or both, then B contains the second part of the complex
          Schur form of the "balanced" versions of the input A and B.
[in]LDB
          LDB is INTEGER
          The leading dimension of B.  LDB >= max(1,N).
[out]ALPHA
          ALPHA is COMPLEX*16 array, dimension (N)
[out]BETA
          BETA is COMPLEX*16 array, dimension (N)
          On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the generalized
          eigenvalues.

          Note: the quotient ALPHA(j)/BETA(j) ) may easily over- or
          underflow, and BETA(j) may even be zero.  Thus, the user
          should avoid naively computing the ratio ALPHA/BETA.
          However, ALPHA 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]VL
          VL is COMPLEX*16 array, dimension (LDVL,N)
          If JOBVL = 'V', the left generalized eigenvectors u(j) are
          stored one after another in the columns of VL, in the same
          order as their eigenvalues.
          Each eigenvector will be scaled so the largest component
          will have abs(real part) + abs(imag. part) = 1.
          Not referenced if JOBVL = 'N'.
[in]LDVL
          LDVL is INTEGER
          The leading dimension of the matrix VL. LDVL >= 1, and
          if JOBVL = 'V', LDVL >= N.
[out]VR
          VR is COMPLEX*16 array, dimension (LDVR,N)
          If JOBVR = 'V', the right generalized eigenvectors v(j) are
          stored one after another in the columns of VR, in the same
          order as their eigenvalues.
          Each eigenvector will be scaled so the largest component
          will have abs(real part) + abs(imag. part) = 1.
          Not referenced if JOBVR = 'N'.
[in]LDVR
          LDVR is INTEGER
          The leading dimension of the matrix VR. LDVR >= 1, and
          if JOBVR = 'V', LDVR >= N.
[out]ILO
          ILO is INTEGER
[out]IHI
          IHI is INTEGER
          ILO and IHI are integer values such that on exit
          A(i,j) = 0 and B(i,j) = 0 if i > j and
          j = 1,...,ILO-1 or i = IHI+1,...,N.
          If BALANC = 'N' or 'S', ILO = 1 and IHI = N.
[out]LSCALE
          LSCALE is DOUBLE PRECISION array, dimension (N)
          Details of the permutations and scaling factors applied
          to the left side of A and B.  If PL(j) is the index of the
          row interchanged with row j, and DL(j) is the scaling
          factor applied to row j, then
            LSCALE(j) = PL(j)  for j = 1,...,ILO-1
                      = DL(j)  for j = ILO,...,IHI
                      = PL(j)  for j = IHI+1,...,N.
          The order in which the interchanges are made is N to IHI+1,
          then 1 to ILO-1.
[out]RSCALE
          RSCALE is DOUBLE PRECISION array, dimension (N)
          Details of the permutations and scaling factors applied
          to the right side of A and B.  If PR(j) is the index of the
          column interchanged with column j, and DR(j) is the scaling
          factor applied to column j, then
            RSCALE(j) = PR(j)  for j = 1,...,ILO-1
                      = DR(j)  for j = ILO,...,IHI
                      = PR(j)  for j = IHI+1,...,N
          The order in which the interchanges are made is N to IHI+1,
          then 1 to ILO-1.
[out]ABNRM
          ABNRM is DOUBLE PRECISION
          The one-norm of the balanced matrix A.
[out]BBNRM
          BBNRM is DOUBLE PRECISION
          The one-norm of the balanced matrix B.
[out]RCONDE
          RCONDE is DOUBLE PRECISION array, dimension (N)
          If SENSE = 'E' or 'B', the reciprocal condition numbers of
          the eigenvalues, stored in consecutive elements of the array.
          If SENSE = 'N' or 'V', RCONDE is not referenced.
[out]RCONDV
          RCONDV is DOUBLE PRECISION array, dimension (N)
          If JOB = 'V' or 'B', the estimated reciprocal condition
          numbers of the eigenvectors, stored in consecutive elements
          of the array. If the eigenvalues cannot be reordered to
          compute RCONDV(j), RCONDV(j) is set to 0; this can only occur
          when the true value would be very small anyway.
          If SENSE = 'N' or 'E', RCONDV is not referenced.
[out]WORK
          WORK is COMPLEX*16 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. LWORK >= max(1,2*N).
          If SENSE = 'E', LWORK >= max(1,4*N).
          If SENSE = 'V' or 'B', LWORK >= max(1,2*N*N+2*N).

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (lrwork)
          lrwork must be at least max(1,6*N) if BALANC = 'S' or 'B',
          and at least max(1,2*N) otherwise.
          Real workspace.
[out]IWORK
          IWORK is INTEGER array, dimension (N+2)
          If SENSE = 'E', IWORK is not referenced.
[out]BWORK
          BWORK is LOGICAL array, dimension (N)
          If SENSE = 'N', BWORK is not referenced.
[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.  No eigenvectors have been
                calculated, but ALPHA(j) and BETA(j) should be correct
                for j=INFO+1,...,N.
          > N:  =N+1: other than QZ iteration failed in ZHGEQZ.
                =N+2: error return from ZTGEVC.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  Balancing a matrix pair (A,B) includes, first, permuting rows and
  columns to isolate eigenvalues, second, applying diagonal similarity
  transformation to the rows and columns to make the rows and columns
  as close in norm as possible. The computed reciprocal condition
  numbers correspond to the balanced matrix. Permuting rows and columns
  will not change the condition numbers (in exact arithmetic) but
  diagonal scaling will.  For further explanation of balancing, see
  section 4.11.1.2 of LAPACK Users' Guide.

  An approximate error bound on the chordal distance between the i-th
  computed generalized eigenvalue w and the corresponding exact
  eigenvalue lambda is

       chord(w, lambda) <= EPS * norm(ABNRM, BBNRM) / RCONDE(I)

  An approximate error bound for the angle between the i-th computed
  eigenvector VL(i) or VR(i) is given by

       EPS * norm(ABNRM, BBNRM) / DIF(i).

  For further explanation of the reciprocal condition numbers RCONDE
  and RCONDV, see section 4.11 of LAPACK User's Guide.

Definition at line 370 of file zggevx.f.

374*
375* -- LAPACK driver routine --
376* -- LAPACK is a software package provided by Univ. of Tennessee, --
377* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
378*
379* .. Scalar Arguments ..
380 CHARACTER BALANC, JOBVL, JOBVR, SENSE
381 INTEGER IHI, ILO, INFO, LDA, LDB, LDVL, LDVR, LWORK, N
382 DOUBLE PRECISION ABNRM, BBNRM
383* ..
384* .. Array Arguments ..
385 LOGICAL BWORK( * )
386 INTEGER IWORK( * )
387 DOUBLE PRECISION LSCALE( * ), RCONDE( * ), RCONDV( * ),
388 $ RSCALE( * ), RWORK( * )
389 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
390 $ BETA( * ), VL( LDVL, * ), VR( LDVR, * ),
391 $ WORK( * )
392* ..
393*
394* =====================================================================
395*
396* .. Parameters ..
397 DOUBLE PRECISION ZERO, ONE
398 parameter( zero = 0.0d+0, one = 1.0d+0 )
399 COMPLEX*16 CZERO, CONE
400 parameter( czero = ( 0.0d+0, 0.0d+0 ),
401 $ cone = ( 1.0d+0, 0.0d+0 ) )
402* ..
403* .. Local Scalars ..
404 LOGICAL ILASCL, ILBSCL, ILV, ILVL, ILVR, LQUERY, NOSCL,
405 $ WANTSB, WANTSE, WANTSN, WANTSV
406 CHARACTER CHTEMP
407 INTEGER I, ICOLS, IERR, IJOBVL, IJOBVR, IN, IROWS,
408 $ ITAU, IWRK, IWRK1, J, JC, JR, M, MAXWRK, MINWRK
409 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
410 $ SMLNUM, TEMP
411 COMPLEX*16 X
412* ..
413* .. Local Arrays ..
414 LOGICAL LDUMMA( 1 )
415* ..
416* .. External Subroutines ..
417 EXTERNAL dlascl, xerbla, zgeqrf, zggbak, zggbal, zgghrd,
419 $ zungqr, zunmqr
420* ..
421* .. External Functions ..
422 LOGICAL LSAME
423 INTEGER ILAENV
424 DOUBLE PRECISION DLAMCH, ZLANGE
425 EXTERNAL lsame, ilaenv, dlamch, zlange
426* ..
427* .. Intrinsic Functions ..
428 INTRINSIC abs, dble, dimag, max, sqrt
429* ..
430* .. Statement Functions ..
431 DOUBLE PRECISION ABS1
432* ..
433* .. Statement Function definitions ..
434 abs1( x ) = abs( dble( x ) ) + abs( dimag( x ) )
435* ..
436* .. Executable Statements ..
437*
438* Decode the input arguments
439*
440 IF( lsame( jobvl, 'N' ) ) THEN
441 ijobvl = 1
442 ilvl = .false.
443 ELSE IF( lsame( jobvl, 'V' ) ) THEN
444 ijobvl = 2
445 ilvl = .true.
446 ELSE
447 ijobvl = -1
448 ilvl = .false.
449 END IF
450*
451 IF( lsame( jobvr, 'N' ) ) THEN
452 ijobvr = 1
453 ilvr = .false.
454 ELSE IF( lsame( jobvr, 'V' ) ) THEN
455 ijobvr = 2
456 ilvr = .true.
457 ELSE
458 ijobvr = -1
459 ilvr = .false.
460 END IF
461 ilv = ilvl .OR. ilvr
462*
463 noscl = lsame( balanc, 'N' ) .OR. lsame( balanc, 'P' )
464 wantsn = lsame( sense, 'N' )
465 wantse = lsame( sense, 'E' )
466 wantsv = lsame( sense, 'V' )
467 wantsb = lsame( sense, 'B' )
468*
469* Test the input arguments
470*
471 info = 0
472 lquery = ( lwork.EQ.-1 )
473 IF( .NOT.( noscl .OR. lsame( balanc,'S' ) .OR.
474 $ lsame( balanc, 'B' ) ) ) THEN
475 info = -1
476 ELSE IF( ijobvl.LE.0 ) THEN
477 info = -2
478 ELSE IF( ijobvr.LE.0 ) THEN
479 info = -3
480 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsb .OR. wantsv ) )
481 $ THEN
482 info = -4
483 ELSE IF( n.LT.0 ) THEN
484 info = -5
485 ELSE IF( lda.LT.max( 1, n ) ) THEN
486 info = -7
487 ELSE IF( ldb.LT.max( 1, n ) ) THEN
488 info = -9
489 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
490 info = -13
491 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
492 info = -15
493 END IF
494*
495* Compute workspace
496* (Note: Comments in the code beginning "Workspace:" describe the
497* minimal amount of workspace needed at that point in the code,
498* as well as the preferred amount for good performance.
499* NB refers to the optimal block size for the immediately
500* following subroutine, as returned by ILAENV. The workspace is
501* computed assuming ILO = 1 and IHI = N, the worst case.)
502*
503 IF( info.EQ.0 ) THEN
504 IF( n.EQ.0 ) THEN
505 minwrk = 1
506 maxwrk = 1
507 ELSE
508 minwrk = 2*n
509 IF( wantse ) THEN
510 minwrk = 4*n
511 ELSE IF( wantsv .OR. wantsb ) THEN
512 minwrk = 2*n*( n + 1)
513 END IF
514 maxwrk = minwrk
515 maxwrk = max( maxwrk,
516 $ n + n*ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
517 maxwrk = max( maxwrk,
518 $ n + n*ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, 0 ) )
519 IF( ilvl ) THEN
520 maxwrk = max( maxwrk, n +
521 $ n*ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, 0 ) )
522 END IF
523 END IF
524 work( 1 ) = maxwrk
525*
526 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
527 info = -25
528 END IF
529 END IF
530*
531 IF( info.NE.0 ) THEN
532 CALL xerbla( 'ZGGEVX', -info )
533 RETURN
534 ELSE IF( lquery ) THEN
535 RETURN
536 END IF
537*
538* Quick return if possible
539*
540 IF( n.EQ.0 )
541 $ RETURN
542*
543* Get machine constants
544*
545 eps = dlamch( 'P' )
546 smlnum = dlamch( 'S' )
547 bignum = one / smlnum
548 smlnum = sqrt( smlnum ) / eps
549 bignum = one / smlnum
550*
551* Scale A if max element outside range [SMLNUM,BIGNUM]
552*
553 anrm = zlange( 'M', n, n, a, lda, rwork )
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 zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
564*
565* Scale B if max element outside range [SMLNUM,BIGNUM]
566*
567 bnrm = zlange( 'M', n, n, b, ldb, rwork )
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 zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
578*
579* Permute and/or balance the matrix pair (A,B)
580* (Real Workspace: need 6*N if BALANC = 'S' or 'B', 1 otherwise)
581*
582 CALL zggbal( balanc, n, a, lda, b, ldb, ilo, ihi, lscale, rscale,
583 $ rwork, ierr )
584*
585* Compute ABNRM and BBNRM
586*
587 abnrm = zlange( '1', n, n, a, lda, rwork( 1 ) )
588 IF( ilascl ) THEN
589 rwork( 1 ) = abnrm
590 CALL dlascl( 'G', 0, 0, anrmto, anrm, 1, 1, rwork( 1 ), 1,
591 $ ierr )
592 abnrm = rwork( 1 )
593 END IF
594*
595 bbnrm = zlange( '1', n, n, b, ldb, rwork( 1 ) )
596 IF( ilbscl ) THEN
597 rwork( 1 ) = bbnrm
598 CALL dlascl( 'G', 0, 0, bnrmto, bnrm, 1, 1, rwork( 1 ), 1,
599 $ ierr )
600 bbnrm = rwork( 1 )
601 END IF
602*
603* Reduce B to triangular form (QR decomposition of B)
604* (Complex Workspace: need N, prefer N*NB )
605*
606 irows = ihi + 1 - ilo
607 IF( ilv .OR. .NOT.wantsn ) THEN
608 icols = n + 1 - ilo
609 ELSE
610 icols = irows
611 END IF
612 itau = 1
613 iwrk = itau + irows
614 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
615 $ work( iwrk ), lwork+1-iwrk, ierr )
616*
617* Apply the unitary transformation to A
618* (Complex Workspace: need N, prefer N*NB)
619*
620 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
621 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
622 $ lwork+1-iwrk, ierr )
623*
624* Initialize VL and/or VR
625* (Workspace: need N, prefer N*NB)
626*
627 IF( ilvl ) THEN
628 CALL zlaset( 'Full', n, n, czero, cone, vl, ldvl )
629 IF( irows.GT.1 ) THEN
630 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
631 $ vl( ilo+1, ilo ), ldvl )
632 END IF
633 CALL zungqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
634 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
635 END IF
636*
637 IF( ilvr )
638 $ CALL zlaset( 'Full', n, n, czero, cone, vr, ldvr )
639*
640* Reduce to generalized Hessenberg form
641* (Workspace: none needed)
642*
643 IF( ilv .OR. .NOT.wantsn ) THEN
644*
645* Eigenvectors requested -- work on whole matrix.
646*
647 CALL zgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
648 $ ldvl, vr, ldvr, ierr )
649 ELSE
650 CALL zgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
651 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
652 END IF
653*
654* Perform QZ algorithm (Compute eigenvalues, and optionally, the
655* Schur forms and Schur vectors)
656* (Complex Workspace: need N)
657* (Real Workspace: need N)
658*
659 iwrk = itau
660 IF( ilv .OR. .NOT.wantsn ) THEN
661 chtemp = 'S'
662 ELSE
663 chtemp = 'E'
664 END IF
665*
666 CALL zhgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
667 $ alpha, beta, vl, ldvl, vr, ldvr, work( iwrk ),
668 $ lwork+1-iwrk, rwork, ierr )
669 IF( ierr.NE.0 ) THEN
670 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
671 info = ierr
672 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
673 info = ierr - n
674 ELSE
675 info = n + 1
676 END IF
677 GO TO 90
678 END IF
679*
680* Compute Eigenvectors and estimate condition numbers if desired
681* ZTGEVC: (Complex Workspace: need 2*N )
682* (Real Workspace: need 2*N )
683* ZTGSNA: (Complex Workspace: need 2*N*N if SENSE='V' or 'B')
684* (Integer Workspace: need N+2 )
685*
686 IF( ilv .OR. .NOT.wantsn ) THEN
687 IF( ilv ) THEN
688 IF( ilvl ) THEN
689 IF( ilvr ) THEN
690 chtemp = 'B'
691 ELSE
692 chtemp = 'L'
693 END IF
694 ELSE
695 chtemp = 'R'
696 END IF
697*
698 CALL ztgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl,
699 $ ldvl, vr, ldvr, n, in, work( iwrk ), rwork,
700 $ ierr )
701 IF( ierr.NE.0 ) THEN
702 info = n + 2
703 GO TO 90
704 END IF
705 END IF
706*
707 IF( .NOT.wantsn ) THEN
708*
709* compute eigenvectors (ZTGEVC) and estimate condition
710* numbers (ZTGSNA). Note that the definition of the condition
711* number is not invariant under transformation (u,v) to
712* (Q*u, Z*v), where (u,v) are eigenvectors of the generalized
713* Schur form (S,T), Q and Z are orthogonal matrices. In order
714* to avoid using extra 2*N*N workspace, we have to
715* re-calculate eigenvectors and estimate the condition numbers
716* one at a time.
717*
718 DO 20 i = 1, n
719*
720 DO 10 j = 1, n
721 bwork( j ) = .false.
722 10 CONTINUE
723 bwork( i ) = .true.
724*
725 iwrk = n + 1
726 iwrk1 = iwrk + n
727*
728 IF( wantse .OR. wantsb ) THEN
729 CALL ztgevc( 'B', 'S', bwork, n, a, lda, b, ldb,
730 $ work( 1 ), n, work( iwrk ), n, 1, m,
731 $ work( iwrk1 ), rwork, ierr )
732 IF( ierr.NE.0 ) THEN
733 info = n + 2
734 GO TO 90
735 END IF
736 END IF
737*
738 CALL ztgsna( sense, 'S', bwork, n, a, lda, b, ldb,
739 $ work( 1 ), n, work( iwrk ), n, rconde( i ),
740 $ rcondv( i ), 1, m, work( iwrk1 ),
741 $ lwork-iwrk1+1, iwork, ierr )
742*
743 20 CONTINUE
744 END IF
745 END IF
746*
747* Undo balancing on VL and VR and normalization
748* (Workspace: none needed)
749*
750 IF( ilvl ) THEN
751 CALL zggbak( balanc, 'L', n, ilo, ihi, lscale, rscale, n, vl,
752 $ ldvl, ierr )
753*
754 DO 50 jc = 1, n
755 temp = zero
756 DO 30 jr = 1, n
757 temp = max( temp, abs1( vl( jr, jc ) ) )
758 30 CONTINUE
759 IF( temp.LT.smlnum )
760 $ GO TO 50
761 temp = one / temp
762 DO 40 jr = 1, n
763 vl( jr, jc ) = vl( jr, jc )*temp
764 40 CONTINUE
765 50 CONTINUE
766 END IF
767*
768 IF( ilvr ) THEN
769 CALL zggbak( balanc, 'R', n, ilo, ihi, lscale, rscale, n, vr,
770 $ ldvr, ierr )
771 DO 80 jc = 1, n
772 temp = zero
773 DO 60 jr = 1, n
774 temp = max( temp, abs1( vr( jr, jc ) ) )
775 60 CONTINUE
776 IF( temp.LT.smlnum )
777 $ GO TO 80
778 temp = one / temp
779 DO 70 jr = 1, n
780 vr( jr, jc ) = vr( jr, jc )*temp
781 70 CONTINUE
782 80 CONTINUE
783 END IF
784*
785* Undo scaling if necessary
786*
787 90 CONTINUE
788*
789 IF( ilascl )
790 $ CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
791*
792 IF( ilbscl )
793 $ CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
794*
795 work( 1 ) = maxwrk
796 RETURN
797*
798* End of ZGGEVX
799*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqrf(m, n, a, lda, tau, work, lwork, info)
ZGEQRF
Definition zgeqrf.f:146
subroutine zggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
ZGGBAK
Definition zggbak.f:148
subroutine zggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
ZGGBAL
Definition zggbal.f:177
subroutine zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:204
subroutine zhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
ZHGEQZ
Definition zhgeqz.f:284
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:143
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.
Definition dlascl.f:143
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ztgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
ZTGEVC
Definition ztgevc.f:219
subroutine ztgsna(job, howmny, select, n, a, lda, b, ldb, vl, ldvl, vr, ldvr, s, dif, mm, m, work, lwork, iwork, info)
ZTGSNA
Definition ztgsna.f:311
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:128
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:167
Here is the call graph for this function:
Here is the caller graph for this function: