LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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
June 2016
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.6.1) --
455 * -- LAPACK is a software package provided by Univ. of Tennessee, --
456 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
457 * June 2016
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  IF( .NOT.lquery .OR. ijob.NE.0 ) THEN
545  DO 10 k = 1, n
546  IF( pair ) THEN
547  pair = .false.
548  ELSE
549  IF( k.LT.n ) THEN
550  IF( a( k+1, k ).EQ.zero ) THEN
551  IF( SELECT( k ) )
552  $ m = m + 1
553  ELSE
554  pair = .true.
555  IF( SELECT( k ) .OR. SELECT( k+1 ) )
556  $ m = m + 2
557  END IF
558  ELSE
559  IF( SELECT( n ) )
560  $ m = m + 1
561  END IF
562  END IF
563  10 CONTINUE
564  END IF
565 *
566  IF( ijob.EQ.1 .OR. ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
567  lwmin = max( 1, 4*n+16, 2*m*(n-m) )
568  liwmin = max( 1, n+6 )
569  ELSE IF( ijob.EQ.3 .OR. ijob.EQ.5 ) THEN
570  lwmin = max( 1, 4*n+16, 4*m*(n-m) )
571  liwmin = max( 1, 2*m*(n-m), n+6 )
572  ELSE
573  lwmin = max( 1, 4*n+16 )
574  liwmin = 1
575  END IF
576 *
577  work( 1 ) = lwmin
578  iwork( 1 ) = liwmin
579 *
580  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
581  info = -22
582  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
583  info = -24
584  END IF
585 *
586  IF( info.NE.0 ) THEN
587  CALL xerbla( 'STGSEN', -info )
588  RETURN
589  ELSE IF( lquery ) THEN
590  RETURN
591  END IF
592 *
593 * Quick return if possible.
594 *
595  IF( m.EQ.n .OR. m.EQ.0 ) THEN
596  IF( wantp ) THEN
597  pl = one
598  pr = one
599  END IF
600  IF( wantd ) THEN
601  dscale = zero
602  dsum = one
603  DO 20 i = 1, n
604  CALL slassq( n, a( 1, i ), 1, dscale, dsum )
605  CALL slassq( n, b( 1, i ), 1, dscale, dsum )
606  20 CONTINUE
607  dif( 1 ) = dscale*sqrt( dsum )
608  dif( 2 ) = dif( 1 )
609  END IF
610  GO TO 60
611  END IF
612 *
613 * Collect the selected blocks at the top-left corner of (A, B).
614 *
615  ks = 0
616  pair = .false.
617  DO 30 k = 1, n
618  IF( pair ) THEN
619  pair = .false.
620  ELSE
621 *
622  swap = SELECT( k )
623  IF( k.LT.n ) THEN
624  IF( a( k+1, k ).NE.zero ) THEN
625  pair = .true.
626  swap = swap .OR. SELECT( k+1 )
627  END IF
628  END IF
629 *
630  IF( swap ) THEN
631  ks = ks + 1
632 *
633 * Swap the K-th block to position KS.
634 * Perform the reordering of diagonal blocks in (A, B)
635 * by orthogonal transformation matrices and update
636 * Q and Z accordingly (if requested):
637 *
638  kk = k
639  IF( k.NE.ks )
640  $ CALL stgexc( wantq, wantz, n, a, lda, b, ldb, q, ldq,
641  $ z, ldz, kk, ks, work, lwork, ierr )
642 *
643  IF( ierr.GT.0 ) THEN
644 *
645 * Swap is rejected: exit.
646 *
647  info = 1
648  IF( wantp ) THEN
649  pl = zero
650  pr = zero
651  END IF
652  IF( wantd ) THEN
653  dif( 1 ) = zero
654  dif( 2 ) = zero
655  END IF
656  GO TO 60
657  END IF
658 *
659  IF( pair )
660  $ ks = ks + 1
661  END IF
662  END IF
663  30 CONTINUE
664  IF( wantp ) THEN
665 *
666 * Solve generalized Sylvester equation for R and L
667 * and compute PL and PR.
668 *
669  n1 = m
670  n2 = n - m
671  i = n1 + 1
672  ijb = 0
673  CALL slacpy( 'Full', n1, n2, a( 1, i ), lda, work, n1 )
674  CALL slacpy( 'Full', n1, n2, b( 1, i ), ldb, work( n1*n2+1 ),
675  $ n1 )
676  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
677  $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ), n1,
678  $ dscale, dif( 1 ), work( n1*n2*2+1 ),
679  $ lwork-2*n1*n2, iwork, ierr )
680 *
681 * Estimate the reciprocal of norms of "projections" onto left
682 * and right eigenspaces.
683 *
684  rdscal = zero
685  dsum = one
686  CALL slassq( n1*n2, work, 1, rdscal, dsum )
687  pl = rdscal*sqrt( dsum )
688  IF( pl.EQ.zero ) THEN
689  pl = one
690  ELSE
691  pl = dscale / ( sqrt( dscale*dscale / pl+pl )*sqrt( pl ) )
692  END IF
693  rdscal = zero
694  dsum = one
695  CALL slassq( n1*n2, work( n1*n2+1 ), 1, rdscal, dsum )
696  pr = rdscal*sqrt( dsum )
697  IF( pr.EQ.zero ) THEN
698  pr = one
699  ELSE
700  pr = dscale / ( sqrt( dscale*dscale / pr+pr )*sqrt( pr ) )
701  END IF
702  END IF
703 *
704  IF( wantd ) THEN
705 *
706 * Compute estimates of Difu and Difl.
707 *
708  IF( wantd1 ) THEN
709  n1 = m
710  n2 = n - m
711  i = n1 + 1
712  ijb = idifjb
713 *
714 * Frobenius norm-based Difu-estimate.
715 *
716  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
717  $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ),
718  $ n1, dscale, dif( 1 ), work( 2*n1*n2+1 ),
719  $ lwork-2*n1*n2, iwork, ierr )
720 *
721 * Frobenius norm-based Difl-estimate.
722 *
723  CALL stgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda, work,
724  $ n2, b( i, i ), ldb, b, ldb, work( n1*n2+1 ),
725  $ n2, dscale, dif( 2 ), work( 2*n1*n2+1 ),
726  $ lwork-2*n1*n2, iwork, ierr )
727  ELSE
728 *
729 *
730 * Compute 1-norm-based estimates of Difu and Difl using
731 * reversed communication with SLACN2. In each step a
732 * generalized Sylvester equation or a transposed variant
733 * is solved.
734 *
735  kase = 0
736  n1 = m
737  n2 = n - m
738  i = n1 + 1
739  ijb = 0
740  mn2 = 2*n1*n2
741 *
742 * 1-norm-based estimate of Difu.
743 *
744  40 CONTINUE
745  CALL slacn2( mn2, work( mn2+1 ), work, iwork, dif( 1 ),
746  $ kase, isave )
747  IF( kase.NE.0 ) THEN
748  IF( kase.EQ.1 ) THEN
749 *
750 * Solve generalized Sylvester equation.
751 *
752  CALL stgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda,
753  $ work, n1, b, ldb, b( i, i ), ldb,
754  $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
755  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
756  $ ierr )
757  ELSE
758 *
759 * Solve the transposed variant.
760 *
761  CALL stgsyl( 'T', ijb, n1, n2, a, lda, a( i, i ), lda,
762  $ work, n1, b, ldb, b( i, i ), ldb,
763  $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
764  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
765  $ ierr )
766  END IF
767  GO TO 40
768  END IF
769  dif( 1 ) = dscale / dif( 1 )
770 *
771 * 1-norm-based estimate of Difl.
772 *
773  50 CONTINUE
774  CALL slacn2( mn2, work( mn2+1 ), work, iwork, dif( 2 ),
775  $ kase, isave )
776  IF( kase.NE.0 ) THEN
777  IF( kase.EQ.1 ) THEN
778 *
779 * Solve generalized Sylvester equation.
780 *
781  CALL stgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda,
782  $ work, n2, b( i, i ), ldb, b, ldb,
783  $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
784  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
785  $ ierr )
786  ELSE
787 *
788 * Solve the transposed variant.
789 *
790  CALL stgsyl( 'T', ijb, n2, n1, a( i, i ), lda, a, lda,
791  $ work, n2, b( i, i ), ldb, b, ldb,
792  $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
793  $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
794  $ ierr )
795  END IF
796  GO TO 50
797  END IF
798  dif( 2 ) = dscale / dif( 2 )
799 *
800  END IF
801  END IF
802 *
803  60 CONTINUE
804 *
805 * Compute generalized eigenvalues of reordered pair (A, B) and
806 * normalize the generalized Schur form.
807 *
808  pair = .false.
809  DO 70 k = 1, n
810  IF( pair ) THEN
811  pair = .false.
812  ELSE
813 *
814  IF( k.LT.n ) THEN
815  IF( a( k+1, k ).NE.zero ) THEN
816  pair = .true.
817  END IF
818  END IF
819 *
820  IF( pair ) THEN
821 *
822 * Compute the eigenvalue(s) at position K.
823 *
824  work( 1 ) = a( k, k )
825  work( 2 ) = a( k+1, k )
826  work( 3 ) = a( k, k+1 )
827  work( 4 ) = a( k+1, k+1 )
828  work( 5 ) = b( k, k )
829  work( 6 ) = b( k+1, k )
830  work( 7 ) = b( k, k+1 )
831  work( 8 ) = b( k+1, k+1 )
832  CALL slag2( work, 2, work( 5 ), 2, smlnum*eps, beta( k ),
833  $ beta( k+1 ), alphar( k ), alphar( k+1 ),
834  $ alphai( k ) )
835  alphai( k+1 ) = -alphai( k )
836 *
837  ELSE
838 *
839  IF( sign( one, b( k, k ) ).LT.zero ) THEN
840 *
841 * If B(K,K) is negative, make it positive
842 *
843  DO 80 i = 1, n
844  a( k, i ) = -a( k, i )
845  b( k, i ) = -b( k, i )
846  IF( wantq ) q( i, k ) = -q( i, k )
847  80 CONTINUE
848  END IF
849 *
850  alphar( k ) = a( k, k )
851  alphai( k ) = zero
852  beta( k ) = b( k, k )
853 *
854  END IF
855  END IF
856  70 CONTINUE
857 *
858  work( 1 ) = lwmin
859  iwork( 1 ) = liwmin
860 *
861  RETURN
862 *
863 * End of STGSEN
864 *
subroutine slassq(N, X, INCX, SCALE, SUMSQ)
SLASSQ updates a sum of squares represented in scaled form.
Definition: slassq.f:105
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
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
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 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 stgexc(WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z, LDZ, IFST, ILST, WORK, LWORK, INFO)
STGEXC
Definition: stgexc.f:222
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69

Here is the call graph for this function:

Here is the caller graph for this function: