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

◆ zggesx()

subroutine zggesx ( character jobvsl,
character jobvsr,
character sort,
external selctg,
character sense,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
integer sdim,
complex*16, dimension( * ) alpha,
complex*16, dimension( * ) beta,
complex*16, dimension( ldvsl, * ) vsl,
integer ldvsl,
complex*16, dimension( ldvsr, * ) vsr,
integer ldvsr,
double precision, dimension( 2 ) rconde,
double precision, dimension( 2 ) rcondv,
complex*16, dimension( * ) work,
integer lwork,
double precision, dimension( * ) rwork,
integer, dimension( * ) iwork,
integer liwork,
logical, dimension( * ) bwork,
integer info )

ZGGESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices

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

Purpose:
!>
!> ZGGESX computes for a pair of N-by-N complex nonsymmetric matrices
!> (A,B), the generalized eigenvalues, the complex Schur form (S,T),
!> and, optionally, the left and/or right matrices of Schur vectors (VSL
!> and VSR).  This gives the generalized Schur factorization
!>
!>      (A,B) = ( (VSL) S (VSR)**H, (VSL) T (VSR)**H )
!>
!> where (VSR)**H is the conjugate-transpose of VSR.
!>
!> Optionally, it also orders the eigenvalues so that a selected cluster
!> of eigenvalues appears in the leading diagonal blocks of the upper
!> triangular matrix S and the upper triangular matrix T; computes
!> a reciprocal condition number for the average of the selected
!> eigenvalues (RCONDE); and computes a reciprocal condition number for
!> the right and left deflating subspaces corresponding to the selected
!> eigenvalues (RCONDV). The leading columns of VSL and VSR then form
!> an orthonormal basis for the corresponding left and right eigenspaces
!> (deflating subspaces).
!>
!> A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
!> or a ratio alpha/beta = w, such that  A - w*B is singular.  It is
!> usually represented as the pair (alpha,beta), as there is a
!> reasonable interpretation for beta=0 or for both being zero.
!>
!> A pair of matrices (S,T) is in generalized complex Schur form if T is
!> upper triangular with non-negative diagonal and S is upper
!> triangular.
!> 
Parameters
[in]JOBVSL
!>          JOBVSL is CHARACTER*1
!>          = 'N':  do not compute the left Schur vectors;
!>          = 'V':  compute the left Schur vectors.
!> 
[in]JOBVSR
!>          JOBVSR is CHARACTER*1
!>          = 'N':  do not compute the right Schur vectors;
!>          = 'V':  compute the right Schur vectors.
!> 
[in]SORT
!>          SORT is CHARACTER*1
!>          Specifies whether or not to order the eigenvalues on the
!>          diagonal of the generalized Schur form.
!>          = 'N':  Eigenvalues are not ordered;
!>          = 'S':  Eigenvalues are ordered (see SELCTG).
!> 
[in]SELCTG
!>          SELCTG is a LOGICAL FUNCTION of two COMPLEX*16 arguments
!>          SELCTG must be declared EXTERNAL in the calling subroutine.
!>          If SORT = 'N', SELCTG is not referenced.
!>          If SORT = 'S', SELCTG is used to select eigenvalues to sort
!>          to the top left of the Schur form.
!>          Note that a selected complex eigenvalue may no longer satisfy
!>          SELCTG(ALPHA(j),BETA(j)) = .TRUE. after ordering, since
!>          ordering may change the value of complex eigenvalues
!>          (especially if the eigenvalue is ill-conditioned), in this
!>          case INFO is set to N+3 see INFO below).
!> 
[in]SENSE
!>          SENSE is CHARACTER*1
!>          Determines which reciprocal condition numbers are computed.
!>          = 'N': None are computed;
!>          = 'E': Computed for average of selected eigenvalues only;
!>          = 'V': Computed for selected deflating subspaces only;
!>          = 'B': Computed for both.
!>          If SENSE = 'E', 'V', or 'B', SORT must equal 'S'.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A, B, VSL, and VSR.  N >= 0.
!> 
[in,out]A
!>          A is COMPLEX*16 array, dimension (LDA, N)
!>          On entry, the first of the pair of matrices.
!>          On exit, A has been overwritten by its generalized Schur
!>          form S.
!> 
[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 second of the pair of matrices.
!>          On exit, B has been overwritten by its generalized Schur
!>          form T.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of B.  LDB >= max(1,N).
!> 
[out]SDIM
!>          SDIM is INTEGER
!>          If SORT = 'N', SDIM = 0.
!>          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
!>          for which SELCTG is true.
!> 
[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.  ALPHA(j) and BETA(j),j=1,...,N  are
!>          the diagonals of the complex Schur form (S,T).  BETA(j) will
!>          be non-negative real.
!>
!>          Note: the quotients 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]VSL
!>          VSL is COMPLEX*16 array, dimension (LDVSL,N)
!>          If JOBVSL = 'V', VSL will contain the left Schur vectors.
!>          Not referenced if JOBVSL = 'N'.
!> 
[in]LDVSL
!>          LDVSL is INTEGER
!>          The leading dimension of the matrix VSL. LDVSL >=1, and
!>          if JOBVSL = 'V', LDVSL >= N.
!> 
[out]VSR
!>          VSR is COMPLEX*16 array, dimension (LDVSR,N)
!>          If JOBVSR = 'V', VSR will contain the right Schur vectors.
!>          Not referenced if JOBVSR = 'N'.
!> 
[in]LDVSR
!>          LDVSR is INTEGER
!>          The leading dimension of the matrix VSR. LDVSR >= 1, and
!>          if JOBVSR = 'V', LDVSR >= N.
!> 
[out]RCONDE
!>          RCONDE is DOUBLE PRECISION array, dimension ( 2 )
!>          If SENSE = 'E' or 'B', RCONDE(1) and RCONDE(2) contain the
!>          reciprocal condition numbers for the average of the selected
!>          eigenvalues.
!>          Not referenced if SENSE = 'N' or 'V'.
!> 
[out]RCONDV
!>          RCONDV is DOUBLE PRECISION array, dimension ( 2 )
!>          If SENSE = 'V' or 'B', RCONDV(1) and RCONDV(2) contain the
!>          reciprocal condition number for the selected deflating
!>          subspaces.
!>          Not referenced if SENSE = 'N' or 'E'.
!> 
[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.
!>          If N = 0, LWORK >= 1, else if SENSE = 'E', 'V', or 'B',
!>          LWORK >= MAX(1,2*N,2*SDIM*(N-SDIM)), else
!>          LWORK >= MAX(1,2*N).  Note that 2*SDIM*(N-SDIM) <= N*N/2.
!>          Note also that an error is only returned if
!>          LWORK < MAX(1,2*N), but if SENSE = 'E' or 'V' or 'B' this may
!>          not be large enough.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the bound on the optimal size of the WORK
!>          array and the minimum size of the IWORK array, returns these
!>          values as the first entries of the WORK and IWORK arrays, and
!>          no error message related to LWORK or LIWORK is issued by
!>          XERBLA.
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension ( 8*N )
!>          Real workspace.
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
!>          On exit, if INFO = 0, IWORK(1) returns the minimum LIWORK.
!> 
[in]LIWORK
!>          LIWORK is INTEGER
!>          The dimension of the array IWORK.
!>          If SENSE = 'N' or N = 0, LIWORK >= 1, otherwise
!>          LIWORK >= N+2.
!>
!>          If LIWORK = -1, then a workspace query is assumed; the
!>          routine only calculates the bound on the optimal size of the
!>          WORK array and the minimum size of the IWORK array, returns
!>          these values as the first entries of the WORK and IWORK
!>          arrays, and no error message related to LWORK or LIWORK is
!>          issued by XERBLA.
!> 
[out]BWORK
!>          BWORK is LOGICAL array, dimension (N)
!>          Not referenced if SORT = 'N'.
!> 
[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.  (A,B) are not in Schur
!>                form, 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: after reordering, roundoff changed values of
!>                      some complex eigenvalues so that leading
!>                      eigenvalues in the Generalized Schur form no
!>                      longer satisfy SELCTG=.TRUE.  This could also
!>                      be caused due to scaling.
!>                =N+3: reordering failed in ZTGSEN.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 324 of file zggesx.f.

329*
330* -- LAPACK driver routine --
331* -- LAPACK is a software package provided by Univ. of Tennessee, --
332* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
333*
334* .. Scalar Arguments ..
335 CHARACTER JOBVSL, JOBVSR, SENSE, SORT
336 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LIWORK, LWORK, N,
337 $ SDIM
338* ..
339* .. Array Arguments ..
340 LOGICAL BWORK( * )
341 INTEGER IWORK( * )
342 DOUBLE PRECISION RCONDE( 2 ), RCONDV( 2 ), RWORK( * )
343 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
344 $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ),
345 $ WORK( * )
346* ..
347* .. Function Arguments ..
348 LOGICAL SELCTG
349 EXTERNAL selctg
350* ..
351*
352* =====================================================================
353*
354* .. Parameters ..
355 DOUBLE PRECISION ZERO, ONE
356 parameter( zero = 0.0d+0, one = 1.0d+0 )
357 COMPLEX*16 CZERO, CONE
358 parameter( czero = ( 0.0d+0, 0.0d+0 ),
359 $ cone = ( 1.0d+0, 0.0d+0 ) )
360* ..
361* .. Local Scalars ..
362 LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
363 $ LQUERY, WANTSB, WANTSE, WANTSN, WANTST, WANTSV
364 INTEGER I, ICOLS, IERR, IHI, IJOB, IJOBVL, IJOBVR,
365 $ ILEFT, ILO, IRIGHT, IROWS, IRWRK, ITAU, IWRK,
366 $ LIWMIN, LWRK, MAXWRK, MINWRK
367 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PL,
368 $ PR, SMLNUM
369* ..
370* .. Local Arrays ..
371 DOUBLE PRECISION DIF( 2 )
372* ..
373* .. External Subroutines ..
374 EXTERNAL xerbla, zgeqrf, zggbak, zggbal, zgghrd,
375 $ zhgeqz,
377* ..
378* .. External Functions ..
379 LOGICAL LSAME
380 INTEGER ILAENV
381 DOUBLE PRECISION DLAMCH, ZLANGE
382 EXTERNAL lsame, ilaenv, dlamch, zlange
383* ..
384* .. Intrinsic Functions ..
385 INTRINSIC max, sqrt
386* ..
387* .. Executable Statements ..
388*
389* Decode the input arguments
390*
391 IF( lsame( jobvsl, 'N' ) ) THEN
392 ijobvl = 1
393 ilvsl = .false.
394 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
395 ijobvl = 2
396 ilvsl = .true.
397 ELSE
398 ijobvl = -1
399 ilvsl = .false.
400 END IF
401*
402 IF( lsame( jobvsr, 'N' ) ) THEN
403 ijobvr = 1
404 ilvsr = .false.
405 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
406 ijobvr = 2
407 ilvsr = .true.
408 ELSE
409 ijobvr = -1
410 ilvsr = .false.
411 END IF
412*
413 wantst = lsame( sort, 'S' )
414 wantsn = lsame( sense, 'N' )
415 wantse = lsame( sense, 'E' )
416 wantsv = lsame( sense, 'V' )
417 wantsb = lsame( sense, 'B' )
418 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
419 IF( wantsn ) THEN
420 ijob = 0
421 ELSE IF( wantse ) THEN
422 ijob = 1
423 ELSE IF( wantsv ) THEN
424 ijob = 2
425 ELSE IF( wantsb ) THEN
426 ijob = 4
427 END IF
428*
429* Test the input arguments
430*
431 info = 0
432 IF( ijobvl.LE.0 ) THEN
433 info = -1
434 ELSE IF( ijobvr.LE.0 ) THEN
435 info = -2
436 ELSE IF( ( .NOT.wantst ) .AND.
437 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
438 info = -3
439 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
440 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
441 info = -5
442 ELSE IF( n.LT.0 ) THEN
443 info = -6
444 ELSE IF( lda.LT.max( 1, n ) ) THEN
445 info = -8
446 ELSE IF( ldb.LT.max( 1, n ) ) THEN
447 info = -10
448 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
449 info = -15
450 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
451 info = -17
452 END IF
453*
454* Compute workspace
455* (Note: Comments in the code beginning "Workspace:" describe the
456* minimal amount of workspace needed at that point in the code,
457* as well as the preferred amount for good performance.
458* NB refers to the optimal block size for the immediately
459* following subroutine, as returned by ILAENV.)
460*
461 IF( info.EQ.0 ) THEN
462 IF( n.GT.0) THEN
463 minwrk = 2*n
464 maxwrk = n*(1 + ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
465 maxwrk = max( maxwrk, n*( 1 +
466 $ ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, -1 ) ) )
467 IF( ilvsl ) THEN
468 maxwrk = max( maxwrk, n*( 1 +
469 $ ilaenv( 1, 'ZUNGQR', ' ', n, 1, n,
470 $ -1 ) ) )
471 END IF
472 lwrk = maxwrk
473 IF( ijob.GE.1 )
474 $ lwrk = max( lwrk, n*n/2 )
475 ELSE
476 minwrk = 1
477 maxwrk = 1
478 lwrk = 1
479 END IF
480 work( 1 ) = lwrk
481 IF( wantsn .OR. n.EQ.0 ) THEN
482 liwmin = 1
483 ELSE
484 liwmin = n + 2
485 END IF
486 iwork( 1 ) = liwmin
487*
488 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
489 info = -21
490 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery) THEN
491 info = -24
492 END IF
493 END IF
494*
495 IF( info.NE.0 ) THEN
496 CALL xerbla( 'ZGGESX', -info )
497 RETURN
498 ELSE IF (lquery) THEN
499 RETURN
500 END IF
501*
502* Quick return if possible
503*
504 IF( n.EQ.0 ) THEN
505 sdim = 0
506 RETURN
507 END IF
508*
509* Get machine constants
510*
511 eps = dlamch( 'P' )
512 smlnum = dlamch( 'S' )
513 bignum = one / smlnum
514 smlnum = sqrt( smlnum ) / eps
515 bignum = one / smlnum
516*
517* Scale A if max element outside range [SMLNUM,BIGNUM]
518*
519 anrm = zlange( 'M', n, n, a, lda, rwork )
520 ilascl = .false.
521 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
522 anrmto = smlnum
523 ilascl = .true.
524 ELSE IF( anrm.GT.bignum ) THEN
525 anrmto = bignum
526 ilascl = .true.
527 END IF
528 IF( ilascl )
529 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
530*
531* Scale B if max element outside range [SMLNUM,BIGNUM]
532*
533 bnrm = zlange( 'M', n, n, b, ldb, rwork )
534 ilbscl = .false.
535 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
536 bnrmto = smlnum
537 ilbscl = .true.
538 ELSE IF( bnrm.GT.bignum ) THEN
539 bnrmto = bignum
540 ilbscl = .true.
541 END IF
542 IF( ilbscl )
543 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
544*
545* Permute the matrix to make it more nearly triangular
546* (Real Workspace: need 6*N)
547*
548 ileft = 1
549 iright = n + 1
550 irwrk = iright + n
551 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
552 $ rwork( iright ), rwork( irwrk ), ierr )
553*
554* Reduce B to triangular form (QR decomposition of B)
555* (Complex Workspace: need N, prefer N*NB)
556*
557 irows = ihi + 1 - ilo
558 icols = n + 1 - ilo
559 itau = 1
560 iwrk = itau + irows
561 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
562 $ work( iwrk ), lwork+1-iwrk, ierr )
563*
564* Apply the unitary transformation to matrix A
565* (Complex Workspace: need N, prefer N*NB)
566*
567 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
568 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
569 $ lwork+1-iwrk, ierr )
570*
571* Initialize VSL
572* (Complex Workspace: need N, prefer N*NB)
573*
574 IF( ilvsl ) THEN
575 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
576 IF( irows.GT.1 ) THEN
577 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
578 $ vsl( ilo+1, ilo ), ldvsl )
579 END IF
580 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
581 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
582 END IF
583*
584* Initialize VSR
585*
586 IF( ilvsr )
587 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
588*
589* Reduce to generalized Hessenberg form
590* (Workspace: none needed)
591*
592 CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
593 $ ldvsl, vsr, ldvsr, ierr )
594*
595 sdim = 0
596*
597* Perform QZ algorithm, computing Schur vectors if desired
598* (Complex Workspace: need N)
599* (Real Workspace: need N)
600*
601 iwrk = itau
602 CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
603 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
604 $ lwork+1-iwrk, rwork( irwrk ), ierr )
605 IF( ierr.NE.0 ) THEN
606 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
607 info = ierr
608 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
609 info = ierr - n
610 ELSE
611 info = n + 1
612 END IF
613 GO TO 40
614 END IF
615*
616* Sort eigenvalues ALPHA/BETA and compute the reciprocal of
617* condition number(s)
618*
619 IF( wantst ) THEN
620*
621* Undo scaling on eigenvalues before SELCTGing
622*
623 IF( ilascl )
624 $ CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n,
625 $ ierr )
626 IF( ilbscl )
627 $ CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n,
628 $ ierr )
629*
630* Select eigenvalues
631*
632 DO 10 i = 1, n
633 bwork( i ) = selctg( alpha( i ), beta( i ) )
634 10 CONTINUE
635*
636* Reorder eigenvalues, transform Generalized Schur vectors, and
637* compute reciprocal condition numbers
638* (Complex Workspace: If IJOB >= 1, need MAX(1, 2*SDIM*(N-SDIM))
639* otherwise, need 1 )
640*
641 CALL ztgsen( ijob, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
642 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, sdim, pl, pr,
643 $ dif, work( iwrk ), lwork-iwrk+1, iwork, liwork,
644 $ ierr )
645*
646 IF( ijob.GE.1 )
647 $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
648 IF( ierr.EQ.-21 ) THEN
649*
650* not enough complex workspace
651*
652 info = -21
653 ELSE
654 IF( ijob.EQ.1 .OR. ijob.EQ.4 ) THEN
655 rconde( 1 ) = pl
656 rconde( 2 ) = pr
657 END IF
658 IF( ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
659 rcondv( 1 ) = dif( 1 )
660 rcondv( 2 ) = dif( 2 )
661 END IF
662 IF( ierr.EQ.1 )
663 $ info = n + 3
664 END IF
665*
666 END IF
667*
668* Apply permutation to VSL and VSR
669* (Workspace: none needed)
670*
671 IF( ilvsl )
672 $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
673 $ rwork( iright ), n, vsl, ldvsl, ierr )
674*
675 IF( ilvsr )
676 $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
677 $ rwork( iright ), n, vsr, ldvsr, ierr )
678*
679* Undo scaling
680*
681 IF( ilascl ) THEN
682 CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
683 CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
684 END IF
685*
686 IF( ilbscl ) THEN
687 CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
688 CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
689 END IF
690*
691 IF( wantst ) THEN
692*
693* Check if reordering is correct
694*
695 lastsl = .true.
696 sdim = 0
697 DO 30 i = 1, n
698 cursl = selctg( alpha( i ), beta( i ) )
699 IF( cursl )
700 $ sdim = sdim + 1
701 IF( cursl .AND. .NOT.lastsl )
702 $ info = n + 2
703 lastsl = cursl
704 30 CONTINUE
705*
706 END IF
707*
708 40 CONTINUE
709*
710 work( 1 ) = maxwrk
711 iwork( 1 ) = liwmin
712*
713 RETURN
714*
715* End of ZGGESX
716*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqrf(m, n, a, lda, tau, work, lwork, info)
ZGEQRF
Definition zgeqrf.f:144
subroutine zggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
ZGGBAK
Definition zggbak.f:147
subroutine zggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
ZGGBAL
Definition zggbal.f:175
subroutine zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:203
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:283
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
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:113
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:142
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:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ztgsen(ijob, wantq, wantz, select, n, a, lda, b, ldb, alpha, beta, q, ldq, z, ldz, m, pl, pr, dif, work, lwork, iwork, liwork, info)
ZTGSEN
Definition ztgsen.f:432
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:126
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:165
Here is the call graph for this function:
Here is the caller graph for this function: