LAPACK  3.10.1
LAPACK: Linear Algebra PACKage

◆ 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 326 of file zggesx.f.

330 *
331 * -- LAPACK driver routine --
332 * -- LAPACK is a software package provided by Univ. of Tennessee, --
333 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
334 *
335 * .. Scalar Arguments ..
336  CHARACTER JOBVSL, JOBVSR, SENSE, SORT
337  INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LIWORK, LWORK, N,
338  $ SDIM
339 * ..
340 * .. Array Arguments ..
341  LOGICAL BWORK( * )
342  INTEGER IWORK( * )
343  DOUBLE PRECISION RCONDE( 2 ), RCONDV( 2 ), RWORK( * )
344  COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
345  $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ),
346  $ WORK( * )
347 * ..
348 * .. Function Arguments ..
349  LOGICAL SELCTG
350  EXTERNAL selctg
351 * ..
352 *
353 * =====================================================================
354 *
355 * .. Parameters ..
356  DOUBLE PRECISION ZERO, ONE
357  parameter( zero = 0.0d+0, one = 1.0d+0 )
358  COMPLEX*16 CZERO, CONE
359  parameter( czero = ( 0.0d+0, 0.0d+0 ),
360  $ cone = ( 1.0d+0, 0.0d+0 ) )
361 * ..
362 * .. Local Scalars ..
363  LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
364  $ LQUERY, WANTSB, WANTSE, WANTSN, WANTST, WANTSV
365  INTEGER I, ICOLS, IERR, IHI, IJOB, IJOBVL, IJOBVR,
366  $ ILEFT, ILO, IRIGHT, IROWS, IRWRK, ITAU, IWRK,
367  $ LIWMIN, LWRK, MAXWRK, MINWRK
368  DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PL,
369  $ PR, SMLNUM
370 * ..
371 * .. Local Arrays ..
372  DOUBLE PRECISION DIF( 2 )
373 * ..
374 * .. External Subroutines ..
375  EXTERNAL dlabad, xerbla, zgeqrf, zggbak, zggbal, zgghrd,
377  $ zunmqr
378 * ..
379 * .. External Functions ..
380  LOGICAL LSAME
381  INTEGER ILAENV
382  DOUBLE PRECISION DLAMCH, ZLANGE
383  EXTERNAL lsame, ilaenv, dlamch, zlange
384 * ..
385 * .. Intrinsic Functions ..
386  INTRINSIC max, sqrt
387 * ..
388 * .. Executable Statements ..
389 *
390 * Decode the input arguments
391 *
392  IF( lsame( jobvsl, 'N' ) ) THEN
393  ijobvl = 1
394  ilvsl = .false.
395  ELSE IF( lsame( jobvsl, 'V' ) ) THEN
396  ijobvl = 2
397  ilvsl = .true.
398  ELSE
399  ijobvl = -1
400  ilvsl = .false.
401  END IF
402 *
403  IF( lsame( jobvsr, 'N' ) ) THEN
404  ijobvr = 1
405  ilvsr = .false.
406  ELSE IF( lsame( jobvsr, 'V' ) ) THEN
407  ijobvr = 2
408  ilvsr = .true.
409  ELSE
410  ijobvr = -1
411  ilvsr = .false.
412  END IF
413 *
414  wantst = lsame( sort, 'S' )
415  wantsn = lsame( sense, 'N' )
416  wantse = lsame( sense, 'E' )
417  wantsv = lsame( sense, 'V' )
418  wantsb = lsame( sense, 'B' )
419  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
420  IF( wantsn ) THEN
421  ijob = 0
422  ELSE IF( wantse ) THEN
423  ijob = 1
424  ELSE IF( wantsv ) THEN
425  ijob = 2
426  ELSE IF( wantsb ) THEN
427  ijob = 4
428  END IF
429 *
430 * Test the input arguments
431 *
432  info = 0
433  IF( ijobvl.LE.0 ) THEN
434  info = -1
435  ELSE IF( ijobvr.LE.0 ) THEN
436  info = -2
437  ELSE IF( ( .NOT.wantst ) .AND. ( .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, -1 ) ) )
470  END IF
471  lwrk = maxwrk
472  IF( ijob.GE.1 )
473  $ lwrk = max( lwrk, n*n/2 )
474  ELSE
475  minwrk = 1
476  maxwrk = 1
477  lwrk = 1
478  END IF
479  work( 1 ) = lwrk
480  IF( wantsn .OR. n.EQ.0 ) THEN
481  liwmin = 1
482  ELSE
483  liwmin = n + 2
484  END IF
485  iwork( 1 ) = liwmin
486 *
487  IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
488  info = -21
489  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery) THEN
490  info = -24
491  END IF
492  END IF
493 *
494  IF( info.NE.0 ) THEN
495  CALL xerbla( 'ZGGESX', -info )
496  RETURN
497  ELSE IF (lquery) THEN
498  RETURN
499  END IF
500 *
501 * Quick return if possible
502 *
503  IF( n.EQ.0 ) THEN
504  sdim = 0
505  RETURN
506  END IF
507 *
508 * Get machine constants
509 *
510  eps = dlamch( 'P' )
511  smlnum = dlamch( 'S' )
512  bignum = one / smlnum
513  CALL dlabad( smlnum, bignum )
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, ierr )
625  IF( ilbscl )
626  $ CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
627 *
628 * Select eigenvalues
629 *
630  DO 10 i = 1, n
631  bwork( i ) = selctg( alpha( i ), beta( i ) )
632  10 CONTINUE
633 *
634 * Reorder eigenvalues, transform Generalized Schur vectors, and
635 * compute reciprocal condition numbers
636 * (Complex Workspace: If IJOB >= 1, need MAX(1, 2*SDIM*(N-SDIM))
637 * otherwise, need 1 )
638 *
639  CALL ztgsen( ijob, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
640  $ alpha, beta, vsl, ldvsl, vsr, ldvsr, sdim, pl, pr,
641  $ dif, work( iwrk ), lwork-iwrk+1, iwork, liwork,
642  $ ierr )
643 *
644  IF( ijob.GE.1 )
645  $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
646  IF( ierr.EQ.-21 ) THEN
647 *
648 * not enough complex workspace
649 *
650  info = -21
651  ELSE
652  IF( ijob.EQ.1 .OR. ijob.EQ.4 ) THEN
653  rconde( 1 ) = pl
654  rconde( 2 ) = pr
655  END IF
656  IF( ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
657  rcondv( 1 ) = dif( 1 )
658  rcondv( 2 ) = dif( 2 )
659  END IF
660  IF( ierr.EQ.1 )
661  $ info = n + 3
662  END IF
663 *
664  END IF
665 *
666 * Apply permutation to VSL and VSR
667 * (Workspace: none needed)
668 *
669  IF( ilvsl )
670  $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
671  $ rwork( iright ), n, vsl, ldvsl, ierr )
672 *
673  IF( ilvsr )
674  $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
675  $ rwork( iright ), n, vsr, ldvsr, ierr )
676 *
677 * Undo scaling
678 *
679  IF( ilascl ) THEN
680  CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
681  CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
682  END IF
683 *
684  IF( ilbscl ) THEN
685  CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
686  CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
687  END IF
688 *
689  IF( wantst ) THEN
690 *
691 * Check if reordering is correct
692 *
693  lastsl = .true.
694  sdim = 0
695  DO 30 i = 1, n
696  cursl = selctg( alpha( i ), beta( i ) )
697  IF( cursl )
698  $ sdim = sdim + 1
699  IF( cursl .AND. .NOT.lastsl )
700  $ info = n + 2
701  lastsl = cursl
702  30 CONTINUE
703 *
704  END IF
705 *
706  40 CONTINUE
707 *
708  work( 1 ) = maxwrk
709  iwork( 1 ) = liwmin
710 *
711  RETURN
712 *
713 * End of ZGGESX
714 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine zggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
ZGGBAL
Definition: zggbal.f:177
subroutine zggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
ZGGBAK
Definition: zggbak.f:148
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 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
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 zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
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
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:433
subroutine zungqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGQR
Definition: zungqr.f:128
subroutine zgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
ZGGHRD
Definition: zgghrd.f:204
subroutine zunmqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMQR
Definition: zunmqr.f:167
subroutine zgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGEQRF VARIANT: left-looking Level 3 BLAS of the algorithm.
Definition: zgeqrf.f:151
Here is the call graph for this function:
Here is the caller graph for this function: