LAPACK 3.12.0
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 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 xerbla, zgeqrf, zggbak, zggbal, zgghrd, 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. ( .NOT.lsame( sort, 'N' ) ) ) THEN
437 info = -3
438 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
439 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
440 info = -5
441 ELSE IF( n.LT.0 ) THEN
442 info = -6
443 ELSE IF( lda.LT.max( 1, n ) ) THEN
444 info = -8
445 ELSE IF( ldb.LT.max( 1, n ) ) THEN
446 info = -10
447 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
448 info = -15
449 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
450 info = -17
451 END IF
452*
453* Compute workspace
454* (Note: Comments in the code beginning "Workspace:" describe the
455* minimal amount of workspace needed at that point in the code,
456* as well as the preferred amount for good performance.
457* NB refers to the optimal block size for the immediately
458* following subroutine, as returned by ILAENV.)
459*
460 IF( info.EQ.0 ) THEN
461 IF( n.GT.0) THEN
462 minwrk = 2*n
463 maxwrk = n*(1 + ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
464 maxwrk = max( maxwrk, n*( 1 +
465 $ ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, -1 ) ) )
466 IF( ilvsl ) THEN
467 maxwrk = max( maxwrk, n*( 1 +
468 $ ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, -1 ) ) )
469 END IF
470 lwrk = maxwrk
471 IF( ijob.GE.1 )
472 $ lwrk = max( lwrk, n*n/2 )
473 ELSE
474 minwrk = 1
475 maxwrk = 1
476 lwrk = 1
477 END IF
478 work( 1 ) = lwrk
479 IF( wantsn .OR. n.EQ.0 ) THEN
480 liwmin = 1
481 ELSE
482 liwmin = n + 2
483 END IF
484 iwork( 1 ) = liwmin
485*
486 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
487 info = -21
488 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery) THEN
489 info = -24
490 END IF
491 END IF
492*
493 IF( info.NE.0 ) THEN
494 CALL xerbla( 'ZGGESX', -info )
495 RETURN
496 ELSE IF (lquery) THEN
497 RETURN
498 END IF
499*
500* Quick return if possible
501*
502 IF( n.EQ.0 ) THEN
503 sdim = 0
504 RETURN
505 END IF
506*
507* Get machine constants
508*
509 eps = dlamch( 'P' )
510 smlnum = dlamch( 'S' )
511 bignum = one / smlnum
512 smlnum = sqrt( smlnum ) / eps
513 bignum = one / smlnum
514*
515* Scale A if max element outside range [SMLNUM,BIGNUM]
516*
517 anrm = zlange( 'M', n, n, a, lda, rwork )
518 ilascl = .false.
519 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
520 anrmto = smlnum
521 ilascl = .true.
522 ELSE IF( anrm.GT.bignum ) THEN
523 anrmto = bignum
524 ilascl = .true.
525 END IF
526 IF( ilascl )
527 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
528*
529* Scale B if max element outside range [SMLNUM,BIGNUM]
530*
531 bnrm = zlange( 'M', n, n, b, ldb, rwork )
532 ilbscl = .false.
533 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
534 bnrmto = smlnum
535 ilbscl = .true.
536 ELSE IF( bnrm.GT.bignum ) THEN
537 bnrmto = bignum
538 ilbscl = .true.
539 END IF
540 IF( ilbscl )
541 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
542*
543* Permute the matrix to make it more nearly triangular
544* (Real Workspace: need 6*N)
545*
546 ileft = 1
547 iright = n + 1
548 irwrk = iright + n
549 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
550 $ rwork( iright ), rwork( irwrk ), ierr )
551*
552* Reduce B to triangular form (QR decomposition of B)
553* (Complex Workspace: need N, prefer N*NB)
554*
555 irows = ihi + 1 - ilo
556 icols = n + 1 - ilo
557 itau = 1
558 iwrk = itau + irows
559 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
560 $ work( iwrk ), lwork+1-iwrk, ierr )
561*
562* Apply the unitary transformation to matrix A
563* (Complex Workspace: need N, prefer N*NB)
564*
565 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
566 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
567 $ lwork+1-iwrk, ierr )
568*
569* Initialize VSL
570* (Complex Workspace: need N, prefer N*NB)
571*
572 IF( ilvsl ) THEN
573 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
574 IF( irows.GT.1 ) THEN
575 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
576 $ vsl( ilo+1, ilo ), ldvsl )
577 END IF
578 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
579 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
580 END IF
581*
582* Initialize VSR
583*
584 IF( ilvsr )
585 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
586*
587* Reduce to generalized Hessenberg form
588* (Workspace: none needed)
589*
590 CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
591 $ ldvsl, vsr, ldvsr, ierr )
592*
593 sdim = 0
594*
595* Perform QZ algorithm, computing Schur vectors if desired
596* (Complex Workspace: need N)
597* (Real Workspace: need N)
598*
599 iwrk = itau
600 CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
601 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
602 $ lwork+1-iwrk, rwork( irwrk ), ierr )
603 IF( ierr.NE.0 ) THEN
604 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
605 info = ierr
606 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
607 info = ierr - n
608 ELSE
609 info = n + 1
610 END IF
611 GO TO 40
612 END IF
613*
614* Sort eigenvalues ALPHA/BETA and compute the reciprocal of
615* condition number(s)
616*
617 IF( wantst ) THEN
618*
619* Undo scaling on eigenvalues before SELCTGing
620*
621 IF( ilascl )
622 $ CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
623 IF( ilbscl )
624 $ CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
625*
626* Select eigenvalues
627*
628 DO 10 i = 1, n
629 bwork( i ) = selctg( alpha( i ), beta( i ) )
630 10 CONTINUE
631*
632* Reorder eigenvalues, transform Generalized Schur vectors, and
633* compute reciprocal condition numbers
634* (Complex Workspace: If IJOB >= 1, need MAX(1, 2*SDIM*(N-SDIM))
635* otherwise, need 1 )
636*
637 CALL ztgsen( ijob, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
638 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, sdim, pl, pr,
639 $ dif, work( iwrk ), lwork-iwrk+1, iwork, liwork,
640 $ ierr )
641*
642 IF( ijob.GE.1 )
643 $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
644 IF( ierr.EQ.-21 ) THEN
645*
646* not enough complex workspace
647*
648 info = -21
649 ELSE
650 IF( ijob.EQ.1 .OR. ijob.EQ.4 ) THEN
651 rconde( 1 ) = pl
652 rconde( 2 ) = pr
653 END IF
654 IF( ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
655 rcondv( 1 ) = dif( 1 )
656 rcondv( 2 ) = dif( 2 )
657 END IF
658 IF( ierr.EQ.1 )
659 $ info = n + 3
660 END IF
661*
662 END IF
663*
664* Apply permutation to VSL and VSR
665* (Workspace: none needed)
666*
667 IF( ilvsl )
668 $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
669 $ rwork( iright ), n, vsl, ldvsl, ierr )
670*
671 IF( ilvsr )
672 $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
673 $ rwork( iright ), n, vsr, ldvsr, ierr )
674*
675* Undo scaling
676*
677 IF( ilascl ) THEN
678 CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
679 CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
680 END IF
681*
682 IF( ilbscl ) THEN
683 CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
684 CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
685 END IF
686*
687 IF( wantst ) THEN
688*
689* Check if reordering is correct
690*
691 lastsl = .true.
692 sdim = 0
693 DO 30 i = 1, n
694 cursl = selctg( alpha( i ), beta( i ) )
695 IF( cursl )
696 $ sdim = sdim + 1
697 IF( cursl .AND. .NOT.lastsl )
698 $ info = n + 2
699 lastsl = cursl
700 30 CONTINUE
701*
702 END IF
703*
704 40 CONTINUE
705*
706 work( 1 ) = maxwrk
707 iwork( 1 ) = liwmin
708*
709 RETURN
710*
711* End of ZGGESX
712*
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 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 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 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: