LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dtgsen ( integer  IJOB,
logical  WANTQ,
logical  WANTZ,
logical, dimension( * )  SELECT,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  ALPHAR,
double precision, dimension( * )  ALPHAI,
double precision, dimension( * )  BETA,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
integer  M,
double precision  PL,
double precision  PR,
double precision, dimension( * )  DIF,
double precision, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

DTGSEN

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

Purpose:
 DTGSEN 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 DGGES), i.e. A is block upper
 triangular with 1-by-1 and 2-by-2 diagonal blocks. B is upper
 triangular.

 DTGSEN also computes the generalized eigenvalues

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

 of the reordered matrix pair (A, B).

 Optionally, DTGSEN 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
[out]ALPHAI
          ALPHAI is DOUBLE PRECISION array, dimension (N)
[out]BETA
          BETA is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION
[out]PR
          PR is DOUBLE PRECISION

          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 DOUBLE PRECISION 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 DOUBLE PRECISION 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:
  DTGSEN 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 DLATDF), then the parameter
  IDIFJB (see below) should be changed from 3 to 4 (routine DLATDF
  (IJOB = 2 will be used)). See DTGSYL 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 454 of file dtgsen.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: