LAPACK  3.5.0
LAPACK: Linear Algebra PACKage
stgsen.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine stgsen (IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, M, PL, PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO)
 STGSEN More...
 

Function/Subroutine Documentation

subroutine stgsen ( integer  IJOB,
logical  WANTQ,
logical  WANTZ,
logical, dimension( * )  SELECT,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  ALPHAR,
real, dimension( * )  ALPHAI,
real, dimension( * )  BETA,
real, dimension( ldq, * )  Q,
integer  LDQ,
real, dimension( ldz, * )  Z,
integer  LDZ,
integer  M,
real  PL,
real  PR,
real, dimension( * )  DIF,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

STGSEN

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

Purpose:
 STGSEN reorders the generalized real Schur decomposition of a real
 matrix pair (A, B) (in terms of an orthonormal equivalence trans-
 formation Q**T * (A, B) * Z), so that a selected cluster of eigenvalues
 appears in the leading diagonal blocks of the upper quasi-triangular
 matrix A and the upper triangular B. The leading columns of Q and
 Z form orthonormal bases of the corresponding left and right eigen-
 spaces (deflating subspaces). (A, B) must be in generalized real
 Schur canonical form (as returned by SGGES), i.e. A is block upper
 triangular with 1-by-1 and 2-by-2 diagonal blocks. B is upper
 triangular.

 STGSEN also computes the generalized eigenvalues

             w(j) = (ALPHAR(j) + i*ALPHAI(j))/BETA(j)

 of the reordered matrix pair (A, B).

 Optionally, STGSEN computes the estimates of reciprocal condition
 numbers for eigenvalues and eigenspaces. These are Difu[(A11,B11),
 (A22,B22)] and Difl[(A11,B11), (A22,B22)], i.e. the separation(s)
 between the matrix pairs (A11, B11) and (A22,B22) that correspond to
 the selected cluster and the eigenvalues outside the cluster, resp.,
 and norms of "projections" onto left and right eigenspaces w.r.t.
 the selected cluster in the (1,1)-block.
Parameters
[in]IJOB
          IJOB is INTEGER
          Specifies whether condition numbers are required for the
          cluster of eigenvalues (PL and PR) or the deflating subspaces
          (Difu and Difl):
           =0: Only reorder w.r.t. SELECT. No extras.
           =1: Reciprocal of norms of "projections" onto left and right
               eigenspaces w.r.t. the selected cluster (PL and PR).
           =2: Upper bounds on Difu and Difl. F-norm-based estimate
               (DIF(1:2)).
           =3: Estimate of Difu and Difl. 1-norm-based estimate
               (DIF(1:2)).
               About 5 times as expensive as IJOB = 2.
           =4: Compute PL, PR and DIF (i.e. 0, 1 and 2 above): Economic
               version to get it all.
           =5: Compute PL, PR and DIF (i.e. 0, 1 and 3 above)
[in]WANTQ
          WANTQ is LOGICAL
          .TRUE. : update the left transformation matrix Q;
          .FALSE.: do not update Q.
[in]WANTZ
          WANTZ is LOGICAL
          .TRUE. : update the right transformation matrix Z;
          .FALSE.: do not update Z.
[in]SELECT
          SELECT is LOGICAL array, dimension (N)
          SELECT specifies the eigenvalues in the selected cluster.
          To select a real eigenvalue w(j), SELECT(j) must be set to
          .TRUE.. To select a complex conjugate pair of eigenvalues
          w(j) and w(j+1), corresponding to a 2-by-2 diagonal block,
          either SELECT(j) or SELECT(j+1) or both must be set to
          .TRUE.; a complex conjugate pair of eigenvalues must be
          either both included in the cluster or both excluded.
[in]N
          N is INTEGER
          The order of the matrices A and B. N >= 0.
[in,out]A
          A is REAL array, dimension(LDA,N)
          On entry, the upper quasi-triangular matrix A, with (A, B) in
          generalized real Schur canonical form.
          On exit, A is overwritten by the reordered matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension(LDB,N)
          On entry, the upper triangular matrix B, with (A, B) in
          generalized real Schur canonical form.
          On exit, B is overwritten by the reordered matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).
[out]ALPHAR
          ALPHAR is REAL array, dimension (N)
[out]ALPHAI
          ALPHAI is REAL array, dimension (N)
[out]BETA
          BETA is REAL array, dimension (N)

          On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
          be the generalized eigenvalues.  ALPHAR(j) + ALPHAI(j)*i
          and BETA(j),j=1,...,N  are the diagonals of the complex Schur
          form (S,T) that would result if the 2-by-2 diagonal blocks of
          the real generalized Schur form of (A,B) were further reduced
          to triangular form using complex unitary transformations.
          If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
          positive, then the j-th and (j+1)-st eigenvalues are a
          complex conjugate pair, with ALPHAI(j+1) negative.
[in,out]Q
          Q is REAL array, dimension (LDQ,N)
          On entry, if WANTQ = .TRUE., Q is an N-by-N matrix.
          On exit, Q has been postmultiplied by the left orthogonal
          transformation matrix which reorder (A, B); The leading M
          columns of Q form orthonormal bases for the specified pair of
          left eigenspaces (deflating subspaces).
          If WANTQ = .FALSE., Q is not referenced.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  LDQ >= 1;
          and if WANTQ = .TRUE., LDQ >= N.
[in,out]Z
          Z is REAL array, dimension (LDZ,N)
          On entry, if WANTZ = .TRUE., Z is an N-by-N matrix.
          On exit, Z has been postmultiplied by the left orthogonal
          transformation matrix which reorder (A, B); The leading M
          columns of Z form orthonormal bases for the specified pair of
          left eigenspaces (deflating subspaces).
          If WANTZ = .FALSE., Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z. LDZ >= 1;
          If WANTZ = .TRUE., LDZ >= N.
[out]M
          M is INTEGER
          The dimension of the specified pair of left and right eigen-
          spaces (deflating subspaces). 0 <= M <= N.
[out]PL
          PL is REAL
[out]PR
          PR is REAL

          If IJOB = 1, 4 or 5, PL, PR are lower bounds on the
          reciprocal of the norm of "projections" onto left and right
          eigenspaces with respect to the selected cluster.
          0 < PL, PR <= 1.
          If M = 0 or M = N, PL = PR  = 1.
          If IJOB = 0, 2 or 3, PL and PR are not referenced.
[out]DIF
          DIF is REAL array, dimension (2).
          If IJOB >= 2, DIF(1:2) store the estimates of Difu and Difl.
          If IJOB = 2 or 4, DIF(1:2) are F-norm-based upper bounds on
          Difu and Difl. If IJOB = 3 or 5, DIF(1:2) are 1-norm-based
          estimates of Difu and Difl.
          If M = 0 or N, DIF(1:2) = F-norm([A, B]).
          If IJOB = 0 or 1, DIF is not referenced.
[out]WORK
          WORK is REAL 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 >=  4*N+16.
          If IJOB = 1, 2 or 4, LWORK >= MAX(4*N+16, 2*M*(N-M)).
          If IJOB = 3 or 5, LWORK >= MAX(4*N+16, 4*M*(N-M)).

          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]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK. LIWORK >= 1.
          If IJOB = 1, 2 or 4, LIWORK >=  N+6.
          If IJOB = 3 or 5, LIWORK >= MAX(2*M*(N-M), N+6).

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal size of the IWORK array,
          returns this value as the first entry of the IWORK array, and
          no error message related to LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
            =0: Successful exit.
            <0: If INFO = -i, the i-th argument had an illegal value.
            =1: Reordering of (A, B) failed because the transformed
                matrix pair (A, B) would be too far from generalized
                Schur form; the problem is very ill-conditioned.
                (A, B) may have been partially reordered.
                If requested, 0 is returned in DIF(*), PL and PR.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  STGSEN first collects the selected eigenvalues by computing
  orthogonal U and W that move them to the top left corner of (A, B).
  In other words, the selected eigenvalues are the eigenvalues of
  (A11, B11) in:

              U**T*(A, B)*W = (A11 A12) (B11 B12) n1
                              ( 0  A22),( 0  B22) n2
                                n1  n2    n1  n2

  where N = n1+n2 and U**T means the transpose of U. The first n1 columns
  of U and W span the specified pair of left and right eigenspaces
  (deflating subspaces) of (A, B).

  If (A, B) has been obtained from the generalized real Schur
  decomposition of a matrix pair (C, D) = Q*(A, B)*Z**T, then the
  reordered generalized real Schur form of (C, D) is given by

           (C, D) = (Q*U)*(U**T*(A, B)*W)*(Z*W)**T,

  and the first n1 columns of Q*U and Z*W span the corresponding
  deflating subspaces of (C, D) (Q and Z store Q*U and Z*W, resp.).

  Note that if the selected eigenvalue is sufficiently ill-conditioned,
  then its value may differ significantly from its value before
  reordering.

  The reciprocal condition numbers of the left and right eigenspaces
  spanned by the first n1 columns of U and W (or Q*U and Z*W) may
  be returned in DIF(1:2), corresponding to Difu and Difl, resp.

  The Difu and Difl are defined as:

       Difu[(A11, B11), (A22, B22)] = sigma-min( Zu )
  and
       Difl[(A11, B11), (A22, B22)] = Difu[(A22, B22), (A11, B11)],

  where sigma-min(Zu) is the smallest singular value of the
  (2*n1*n2)-by-(2*n1*n2) matrix

       Zu = [ kron(In2, A11)  -kron(A22**T, In1) ]
            [ kron(In2, B11)  -kron(B22**T, In1) ].

  Here, Inx is the identity matrix of size nx and A22**T is the
  transpose of A22. kron(X, Y) is the Kronecker product between
  the matrices X and Y.

  When DIF(2) is small, small changes in (A, B) can cause large changes
  in the deflating subspace. An approximate (asymptotic) bound on the
  maximum angular error in the computed deflating subspaces is

       EPS * norm((A, B)) / DIF(2),

  where EPS is the machine precision.

  The reciprocal norm of the projectors on the left and right
  eigenspaces associated with (A11, B11) may be returned in PL and PR.
  They are computed as follows. First we compute L and R so that
  P*(A, B)*Q is block diagonal, where

       P = ( I -L ) n1           Q = ( I R ) n1
           ( 0  I ) n2    and        ( 0 I ) n2
             n1 n2                    n1 n2

  and (L, R) is the solution to the generalized Sylvester equation

       A11*R - L*A22 = -A12
       B11*R - L*B22 = -B12

  Then PL = (F-norm(L)**2+1)**(-1/2) and PR = (F-norm(R)**2+1)**(-1/2).
  An approximate (asymptotic) bound on the average absolute error of
  the selected eigenvalues is

       EPS * norm((A, B)) / PL.

  There are also global error bounds which valid for perturbations up
  to a certain restriction:  A lower bound (x) on the smallest
  F-norm(E,F) for which an eigenvalue of (A11, B11) may move and
  coalesce with an eigenvalue of (A22, B22) under perturbation (E,F),
  (i.e. (A + E, B + F), is

   x = min(Difu,Difl)/((1/(PL*PL)+1/(PR*PR))**(1/2)+2*max(1/PL,1/PR)).

  An approximate bound on x can be computed from DIF(1:2), PL and PR.

  If y = ( F-norm(E,F) / x) <= 1, the angles between the perturbed
  (L', R') and unperturbed (L, R) left and right deflating subspaces
  associated with the selected cluster in the (1,1)-blocks can be
  bounded as

   max-angle(L, L') <= arctan( y * PL / (1 - y * (1 - PL * PL)**(1/2))
   max-angle(R, R') <= arctan( y * PR / (1 - y * (1 - PR * PR)**(1/2))

  See LAPACK User's Guide section 4.11 or the following references
  for more information.

  Note that if the default method for computing the Frobenius-norm-
  based estimate DIF is not wanted (see SLATDF), then the parameter
  IDIFJB (see below) should be changed from 3 to 4 (routine SLATDF
  (IJOB = 2 will be used)). See STGSYL for more details.
Contributors:
Bo Kagstrom and Peter Poromaa, Department of Computing Science, Umea University, S-901 87 Umea, Sweden.
References:
  [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
      Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
      M.S. Moonen et al (eds), Linear Algebra for Large Scale and
      Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.

  [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
      Eigenvalues of a Regular Matrix Pair (A, B) and Condition
      Estimation: Theory, Algorithms and Software,
      Report UMINF - 94.04, Department of Computing Science, Umea
      University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
      Note 87. To appear in Numerical Algorithms, 1996.

  [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
      for Solving the Generalized Sylvester Equation and Estimating the
      Separation between Regular Matrix Pairs, Report UMINF - 93.23,
      Department of Computing Science, Umea University, S-901 87 Umea,
      Sweden, December 1993, Revised April 1994, Also as LAPACK Working
      Note 75. To appear in ACM Trans. on Math. Software, Vol 22, No 1,
      1996.

Definition at line 453 of file stgsen.f.

453 *
454 * -- LAPACK computational routine (version 3.4.0) --
455 * -- LAPACK is a software package provided by Univ. of Tennessee, --
456 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
457 * November 2011
458 *
459 * .. Scalar Arguments ..
460  LOGICAL wantq, wantz
461  INTEGER ijob, info, lda, ldb, ldq, ldz, liwork, lwork,
462  $ m, n
463  REAL pl, pr
464 * ..
465 * .. Array Arguments ..
466  LOGICAL select( * )
467  INTEGER iwork( * )
468  REAL a( lda, * ), alphai( * ), alphar( * ),
469  $ b( ldb, * ), beta( * ), dif( * ), q( ldq, * ),
470  $ work( * ), z( ldz, * )
471 * ..
472 *
473 * =====================================================================
474 *
475 * .. Parameters ..
476  INTEGER idifjb
477  parameter( idifjb = 3 )
478  REAL zero, one
479  parameter( zero = 0.0e+0, one = 1.0e+0 )
480 * ..
481 * .. Local Scalars ..
482  LOGICAL lquery, pair, swap, wantd, wantd1, wantd2,
483  $ wantp
484  INTEGER i, ierr, ijb, k, kase, kk, ks, liwmin, lwmin,
485  $ mn2, n1, n2
486  REAL dscale, dsum, eps, rdscal, smlnum
487 * ..
488 * .. Local Arrays ..
489  INTEGER isave( 3 )
490 * ..
491 * .. External Subroutines ..
492  EXTERNAL slacn2, slacpy, slag2, slassq, stgexc, stgsyl,
493  $ xerbla
494 * ..
495 * .. External Functions ..
496  REAL slamch
497  EXTERNAL slamch
498 * ..
499 * .. Intrinsic Functions ..
500  INTRINSIC max, sign, sqrt
501 * ..
502 * .. Executable Statements ..
503 *
504 * Decode and test the input parameters
505 *
506  info = 0
507  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
508 *
509  IF( ijob.LT.0 .OR. ijob.GT.5 ) THEN
510  info = -1
511  ELSE IF( n.LT.0 ) THEN
512  info = -5
513  ELSE IF( lda.LT.max( 1, n ) ) THEN
514  info = -7
515  ELSE IF( ldb.LT.max( 1, n ) ) THEN
516  info = -9
517  ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
518  info = -14
519  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
520  info = -16
521  END IF
522 *
523  IF( info.NE.0 ) THEN
524  CALL xerbla( 'STGSEN', -info )
525  RETURN
526  END IF
527 *
528 * Get machine constants
529 *
530  eps = slamch( 'P' )
531  smlnum = slamch( 'S' ) / eps
532  ierr = 0
533 *
534  wantp = ijob.EQ.1 .OR. ijob.GE.4
535  wantd1 = ijob.EQ.2 .OR. ijob.EQ.4
536  wantd2 = ijob.EQ.3 .OR. ijob.EQ.5
537  wantd = wantd1 .OR. wantd2
538 *
539 * Set M to the dimension of the specified pair of deflating
540 * subspaces.
541 *
542  m = 0
543  pair = .false.
544  DO 10 k = 1, n
545  IF( pair ) THEN
546  pair = .false.
547  ELSE
548  IF( k.LT.n ) THEN
549  IF( a( k+1, k ).EQ.zero ) THEN
550  IF( SELECT( k ) )
551  $ m = m + 1
552  ELSE
553  pair = .true.
554  IF( SELECT( k ) .OR. SELECT( k+1 ) )
555  $ m = m + 2
556  END IF
557  ELSE
558  IF( SELECT( n ) )
559  $ m = m + 1
560  END IF
561  END IF
562  10 CONTINUE
563 *
564  IF( ijob.EQ.1 .OR. ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
565  lwmin = max( 1, 4*n+16, 2*m*(n-m) )
566  liwmin = max( 1, n+6 )
567  ELSE IF( ijob.EQ.3 .OR. ijob.EQ.5 ) THEN
568  lwmin = max( 1, 4*n+16, 4*m*(n-m) )
569  liwmin = max( 1, 2*m*(n-m), n+6 )
570  ELSE
571  lwmin = max( 1, 4*n+16 )
572  liwmin = 1
573  END IF
574 *
575  work( 1 ) = lwmin
576  iwork( 1 ) = liwmin
577 *
578  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
579  info = -22
580  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
581  info = -24
582  END IF
583 *
584  IF( info.NE.0 ) THEN
585  CALL xerbla( 'STGSEN', -info )
586  RETURN
587  ELSE IF( lquery ) THEN
588  RETURN
589  END IF
590 *
591 * Quick return if possible.
592 *
593  IF( m.EQ.n .OR. m.EQ.0 ) THEN
594  IF( wantp ) THEN
595  pl = one
596  pr = one
597  END IF
598  IF( wantd ) THEN
599  dscale = zero
600  dsum = one
601  DO 20 i = 1, n
602  CALL slassq( n, a( 1, i ), 1, dscale, dsum )
603  CALL slassq( n, b( 1, i ), 1, dscale, dsum )
604  20 CONTINUE
605  dif( 1 ) = dscale*sqrt( dsum )
606  dif( 2 ) = dif( 1 )
607  END IF
608  GO TO 60
609  END IF
610 *
611 * Collect the selected blocks at the top-left corner of (A, B).
612 *
613  ks = 0
614  pair = .false.
615  DO 30 k = 1, n
616  IF( pair ) THEN
617  pair = .false.
618  ELSE
619 *
620  swap = SELECT( k )
621  IF( k.LT.n ) THEN
622  IF( a( k+1, k ).NE.zero ) THEN
623  pair = .true.
624  swap = swap .OR. SELECT( k+1 )
625  END IF
626  END IF
627 *
628  IF( swap ) THEN
629  ks = ks + 1
630 *
631 * Swap the K-th block to position KS.
632 * Perform the reordering of diagonal blocks in (A, B)
633 * by orthogonal transformation matrices and update
634 * Q and Z accordingly (if requested):
635 *
636  kk = k
637  IF( k.NE.ks )
638  $ CALL stgexc( wantq, wantz, n, a, lda, b, ldb, q, ldq,
639  $ z, ldz, kk, ks, work, lwork, ierr )
640 *
641  IF( ierr.GT.0 ) THEN
642 *
643 * Swap is rejected: exit.
644 *
645  info = 1
646  IF( wantp ) THEN
647  pl = zero
648  pr = zero
649  END IF
650  IF( wantd ) THEN
651  dif( 1 ) = zero
652  dif( 2 ) = zero
653  END IF
654  GO TO 60
655  END IF
656 *
657  IF( pair )
658  $ ks = ks + 1
659  END IF
660  END IF
661  30 CONTINUE
662  IF( wantp ) THEN
663 *
664 * Solve generalized Sylvester equation for R and L
665 * and compute PL and PR.
666 *
667  n1 = m
668  n2 = n - m
669  i = n1 + 1
670  ijb = 0
671  CALL slacpy( 'Full', n1, n2, a( 1, i ), lda, work, n1 )
672  CALL slacpy( 'Full', n1, n2, b( 1, i ), ldb, work( n1*n2+1 ),
673  $ n1 )
674  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
675  $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ), n1,
676  $ dscale, dif( 1 ), work( n1*n2*2+1 ),
677  $ lwork-2*n1*n2, iwork, ierr )
678 *
679 * Estimate the reciprocal of norms of "projections" onto left
680 * and right eigenspaces.
681 *
682  rdscal = zero
683  dsum = one
684  CALL slassq( n1*n2, work, 1, rdscal, dsum )
685  pl = rdscal*sqrt( dsum )
686  IF( pl.EQ.zero ) THEN
687  pl = one
688  ELSE
689  pl = dscale / ( sqrt( dscale*dscale / pl+pl )*sqrt( pl ) )
690  END IF
691  rdscal = zero
692  dsum = one
693  CALL slassq( n1*n2, work( n1*n2+1 ), 1, rdscal, dsum )
694  pr = rdscal*sqrt( dsum )
695  IF( pr.EQ.zero ) THEN
696  pr = one
697  ELSE
698  pr = dscale / ( sqrt( dscale*dscale / pr+pr )*sqrt( pr ) )
699  END IF
700  END IF
701 *
702  IF( wantd ) THEN
703 *
704 * Compute estimates of Difu and Difl.
705 *
706  IF( wantd1 ) THEN
707  n1 = m
708  n2 = n - m
709  i = n1 + 1
710  ijb = idifjb
711 *
712 * Frobenius norm-based Difu-estimate.
713 *
714  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
715  $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ),
716  $ n1, dscale, dif( 1 ), work( 2*n1*n2+1 ),
717  $ lwork-2*n1*n2, iwork, ierr )
718 *
719 * Frobenius norm-based Difl-estimate.
720 *
721  CALL stgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda, work,
722  $ n2, b( i, i ), ldb, b, ldb, work( n1*n2+1 ),
723  $ n2, dscale, dif( 2 ), work( 2*n1*n2+1 ),
724  $ lwork-2*n1*n2, iwork, ierr )
725  ELSE
726 *
727 *
728 * Compute 1-norm-based estimates of Difu and Difl using
729 * reversed communication with SLACN2. In each step a
730 * generalized Sylvester equation or a transposed variant
731 * is solved.
732 *
733  kase = 0
734  n1 = m
735  n2 = n - m
736  i = n1 + 1
737  ijb = 0
738  mn2 = 2*n1*n2
739 *
740 * 1-norm-based estimate of Difu.
741 *
742  40 CONTINUE
743  CALL slacn2( mn2, work( mn2+1 ), work, iwork, dif( 1 ),
744  $ kase, isave )
745  IF( kase.NE.0 ) THEN
746  IF( kase.EQ.1 ) THEN
747 *
748 * Solve generalized Sylvester equation.
749 *
750  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda,
751  $ work, n1, b, ldb, b( i, i ), ldb,
752  $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
753  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
754  $ ierr )
755  ELSE
756 *
757 * Solve the transposed variant.
758 *
759  CALL stgsyl( 'T', ijb, n1, n2, a, lda, a( i, i ), lda,
760  $ work, n1, b, ldb, b( i, i ), ldb,
761  $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
762  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
763  $ ierr )
764  END IF
765  GO TO 40
766  END IF
767  dif( 1 ) = dscale / dif( 1 )
768 *
769 * 1-norm-based estimate of Difl.
770 *
771  50 CONTINUE
772  CALL slacn2( mn2, work( mn2+1 ), work, iwork, dif( 2 ),
773  $ kase, isave )
774  IF( kase.NE.0 ) THEN
775  IF( kase.EQ.1 ) THEN
776 *
777 * Solve generalized Sylvester equation.
778 *
779  CALL stgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda,
780  $ work, n2, b( i, i ), ldb, b, ldb,
781  $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
782  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
783  $ ierr )
784  ELSE
785 *
786 * Solve the transposed variant.
787 *
788  CALL stgsyl( 'T', ijb, n2, n1, a( i, i ), lda, a, lda,
789  $ work, n2, b( i, i ), ldb, b, ldb,
790  $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
791  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
792  $ ierr )
793  END IF
794  GO TO 50
795  END IF
796  dif( 2 ) = dscale / dif( 2 )
797 *
798  END IF
799  END IF
800 *
801  60 CONTINUE
802 *
803 * Compute generalized eigenvalues of reordered pair (A, B) and
804 * normalize the generalized Schur form.
805 *
806  pair = .false.
807  DO 70 k = 1, n
808  IF( pair ) THEN
809  pair = .false.
810  ELSE
811 *
812  IF( k.LT.n ) THEN
813  IF( a( k+1, k ).NE.zero ) THEN
814  pair = .true.
815  END IF
816  END IF
817 *
818  IF( pair ) THEN
819 *
820 * Compute the eigenvalue(s) at position K.
821 *
822  work( 1 ) = a( k, k )
823  work( 2 ) = a( k+1, k )
824  work( 3 ) = a( k, k+1 )
825  work( 4 ) = a( k+1, k+1 )
826  work( 5 ) = b( k, k )
827  work( 6 ) = b( k+1, k )
828  work( 7 ) = b( k, k+1 )
829  work( 8 ) = b( k+1, k+1 )
830  CALL slag2( work, 2, work( 5 ), 2, smlnum*eps, beta( k ),
831  $ beta( k+1 ), alphar( k ), alphar( k+1 ),
832  $ alphai( k ) )
833  alphai( k+1 ) = -alphai( k )
834 *
835  ELSE
836 *
837  IF( sign( one, b( k, k ) ).LT.zero ) THEN
838 *
839 * If B(K,K) is negative, make it positive
840 *
841  DO 80 i = 1, n
842  a( k, i ) = -a( k, i )
843  b( k, i ) = -b( k, i )
844  IF( wantq ) q( i, k ) = -q( i, k )
845  80 CONTINUE
846  END IF
847 *
848  alphar( k ) = a( k, k )
849  alphai( k ) = zero
850  beta( k ) = b( k, k )
851 *
852  END IF
853  END IF
854  70 CONTINUE
855 *
856  work( 1 ) = lwmin
857  iwork( 1 ) = liwmin
858 *
859  RETURN
860 *
861 * End of STGSEN
862 *
subroutine slassq(N, X, INCX, SCALE, SUMSQ)
SLASSQ updates a sum of squares represented in scaled form.
Definition: slassq.f:105
subroutine slag2(A, LDA, B, LDB, SAFMIN, SCALE1, SCALE2, WR1, WR2, WI)
SLAG2 computes the eigenvalues of a 2-by-2 generalized eigenvalue problem, with scaling as necessary ...
Definition: slag2.f:158
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine slacn2(N, V, X, ISGN, EST, KASE, ISAVE)
SLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: slacn2.f:138
subroutine stgexc(WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z, LDZ, IFST, ILST, WORK, LWORK, INFO)
STGEXC
Definition: stgexc.f:222
subroutine stgsyl(TRANS, IJOB, M, N, A, LDA, B, LDB, C, LDC, D, LDD, E, LDE, F, LDF, SCALE, DIF, WORK, LWORK, IWORK, INFO)
STGSYL
Definition: stgsyl.f:301

Here is the call graph for this function:

Here is the caller graph for this function: