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

◆ cggevx()

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

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

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

Purpose:
!>
!> CGGEVX 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 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  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 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  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 array, dimension (N)
!> 
[out]BETA
!>          BETA is COMPLEX 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 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 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 REAL 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 REAL 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 REAL
!>          The one-norm of the balanced matrix A.
!> 
[out]BBNRM
!>          BBNRM is REAL
!>          The one-norm of the balanced matrix B.
!> 
[out]RCONDE
!>          RCONDE is REAL 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 REAL array, dimension (N)
!>          If SENSE = '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 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 REAL 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 CHGEQZ.
!>                =N+2: error return from CTGEVC.
!> 
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 368 of file cggevx.f.

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