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

◆ cggesx()

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

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

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

Purpose:
!>
!> CGGESX 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 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 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 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 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.  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 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 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 REAL 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 REAL 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 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 REAL 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 WORK.
!>          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 CHGEQZ
!>                =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 CTGSEN.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

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