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

◆ dtgsen()

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  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  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  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.
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 446 of file dtgsen.f.

450*
451* -- LAPACK computational routine --
452* -- LAPACK is a software package provided by Univ. of Tennessee, --
453* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
454*
455* .. Scalar Arguments ..
456 LOGICAL WANTQ, WANTZ
457 INTEGER IJOB, INFO, LDA, LDB, LDQ, LDZ, LIWORK, LWORK,
458 $ M, N
459 DOUBLE PRECISION PL, PR
460* ..
461* .. Array Arguments ..
462 LOGICAL SELECT( * )
463 INTEGER IWORK( * )
464 DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
465 $ B( LDB, * ), BETA( * ), DIF( * ), Q( LDQ, * ),
466 $ WORK( * ), Z( LDZ, * )
467* ..
468*
469* =====================================================================
470*
471* .. Parameters ..
472 INTEGER IDIFJB
473 parameter( idifjb = 3 )
474 DOUBLE PRECISION ZERO, ONE
475 parameter( zero = 0.0d+0, one = 1.0d+0 )
476* ..
477* .. Local Scalars ..
478 LOGICAL LQUERY, PAIR, SWAP, WANTD, WANTD1, WANTD2,
479 $ WANTP
480 INTEGER I, IERR, IJB, K, KASE, KK, KS, LIWMIN, LWMIN,
481 $ MN2, N1, N2
482 DOUBLE PRECISION DSCALE, DSUM, EPS, RDSCAL, SMLNUM
483* ..
484* .. Local Arrays ..
485 INTEGER ISAVE( 3 )
486* ..
487* .. External Subroutines ..
488 EXTERNAL dlacn2, dlacpy, dlag2, dlassq, dtgexc,
489 $ dtgsyl,
490 $ xerbla
491* ..
492* .. External Functions ..
493 DOUBLE PRECISION DLAMCH
494 EXTERNAL dlamch
495* ..
496* .. Intrinsic Functions ..
497 INTRINSIC max, sign, sqrt
498* ..
499* .. Executable Statements ..
500*
501* Decode and test the input parameters
502*
503 info = 0
504 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
505*
506 IF( ijob.LT.0 .OR. ijob.GT.5 ) THEN
507 info = -1
508 ELSE IF( n.LT.0 ) THEN
509 info = -5
510 ELSE IF( lda.LT.max( 1, n ) ) THEN
511 info = -7
512 ELSE IF( ldb.LT.max( 1, n ) ) THEN
513 info = -9
514 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
515 info = -14
516 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
517 info = -16
518 END IF
519*
520 IF( info.NE.0 ) THEN
521 CALL xerbla( 'DTGSEN', -info )
522 RETURN
523 END IF
524*
525* Get machine constants
526*
527 eps = dlamch( 'P' )
528 smlnum = dlamch( 'S' ) / eps
529 ierr = 0
530*
531 wantp = ijob.EQ.1 .OR. ijob.GE.4
532 wantd1 = ijob.EQ.2 .OR. ijob.EQ.4
533 wantd2 = ijob.EQ.3 .OR. ijob.EQ.5
534 wantd = wantd1 .OR. wantd2
535*
536* Set M to the dimension of the specified pair of deflating
537* subspaces.
538*
539 m = 0
540 pair = .false.
541 IF( .NOT.lquery .OR. ijob.NE.0 ) THEN
542 DO 10 k = 1, n
543 IF( pair ) THEN
544 pair = .false.
545 ELSE
546 IF( k.LT.n ) THEN
547 IF( a( k+1, k ).EQ.zero ) THEN
548 IF( SELECT( k ) )
549 $ m = m + 1
550 ELSE
551 pair = .true.
552 IF( SELECT( k ) .OR. SELECT( k+1 ) )
553 $ m = m + 2
554 END IF
555 ELSE
556 IF( SELECT( n ) )
557 $ m = m + 1
558 END IF
559 END IF
560 10 CONTINUE
561 END IF
562*
563 IF( ijob.EQ.1 .OR. ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
564 lwmin = max( 1, 4*n+16, 2*m*( n-m ) )
565 liwmin = max( 1, n+6 )
566 ELSE IF( ijob.EQ.3 .OR. ijob.EQ.5 ) THEN
567 lwmin = max( 1, 4*n+16, 4*m*( n-m ) )
568 liwmin = max( 1, 2*m*( n-m ), n+6 )
569 ELSE
570 lwmin = max( 1, 4*n+16 )
571 liwmin = 1
572 END IF
573*
574 work( 1 ) = lwmin
575 iwork( 1 ) = liwmin
576*
577 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
578 info = -22
579 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
580 info = -24
581 END IF
582*
583 IF( info.NE.0 ) THEN
584 CALL xerbla( 'DTGSEN', -info )
585 RETURN
586 ELSE IF( lquery ) THEN
587 RETURN
588 END IF
589*
590* Quick return if possible.
591*
592 IF( m.EQ.n .OR. m.EQ.0 ) THEN
593 IF( wantp ) THEN
594 pl = one
595 pr = one
596 END IF
597 IF( wantd ) THEN
598 dscale = zero
599 dsum = one
600 DO 20 i = 1, n
601 CALL dlassq( n, a( 1, i ), 1, dscale, dsum )
602 CALL dlassq( n, b( 1, i ), 1, dscale, dsum )
603 20 CONTINUE
604 dif( 1 ) = dscale*sqrt( dsum )
605 dif( 2 ) = dif( 1 )
606 END IF
607 GO TO 60
608 END IF
609*
610* Collect the selected blocks at the top-left corner of (A, B).
611*
612 ks = 0
613 pair = .false.
614 DO 30 k = 1, n
615 IF( pair ) THEN
616 pair = .false.
617 ELSE
618*
619 swap = SELECT( k )
620 IF( k.LT.n ) THEN
621 IF( a( k+1, k ).NE.zero ) THEN
622 pair = .true.
623 swap = swap .OR. SELECT( k+1 )
624 END IF
625 END IF
626*
627 IF( swap ) THEN
628 ks = ks + 1
629*
630* Swap the K-th block to position KS.
631* Perform the reordering of diagonal blocks in (A, B)
632* by orthogonal transformation matrices and update
633* Q and Z accordingly (if requested):
634*
635 kk = k
636 IF( k.NE.ks )
637 $ CALL dtgexc( wantq, wantz, n, a, lda, b, ldb, q,
638 $ 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 dlacpy( 'Full', n1, n2, a( 1, i ), lda, work, n1 )
672 CALL dlacpy( 'Full', n1, n2, b( 1, i ), ldb,
673 $ work( n1*n2+1 ),
674 $ n1 )
675 CALL dtgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
676 $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ), n1,
677 $ dscale, dif( 1 ), work( n1*n2*2+1 ),
678 $ lwork-2*n1*n2, iwork, ierr )
679*
680* Estimate the reciprocal of norms of "projections" onto left
681* and right eigenspaces.
682*
683 rdscal = zero
684 dsum = one
685 CALL dlassq( n1*n2, work, 1, rdscal, dsum )
686 pl = rdscal*sqrt( dsum )
687 IF( pl.EQ.zero ) THEN
688 pl = one
689 ELSE
690 pl = dscale / ( sqrt( dscale*dscale / pl+pl )*sqrt( pl ) )
691 END IF
692 rdscal = zero
693 dsum = one
694 CALL dlassq( n1*n2, work( n1*n2+1 ), 1, rdscal, dsum )
695 pr = rdscal*sqrt( dsum )
696 IF( pr.EQ.zero ) THEN
697 pr = one
698 ELSE
699 pr = dscale / ( sqrt( dscale*dscale / pr+pr )*sqrt( pr ) )
700 END IF
701 END IF
702*
703 IF( wantd ) THEN
704*
705* Compute estimates of Difu and Difl.
706*
707 IF( wantd1 ) THEN
708 n1 = m
709 n2 = n - m
710 i = n1 + 1
711 ijb = idifjb
712*
713* Frobenius norm-based Difu-estimate.
714*
715 CALL dtgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda,
716 $ 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 dtgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda,
724 $ 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 ),
754 $ lda,
755 $ work, n1, b, ldb, b( i, i ), ldb,
756 $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
757 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
758 $ ierr )
759 ELSE
760*
761* Solve the transposed variant.
762*
763 CALL dtgsyl( 'T', ijb, n1, n2, a, lda, a( i, i ),
764 $ lda,
765 $ work, n1, b, ldb, b( i, i ), ldb,
766 $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
767 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
768 $ ierr )
769 END IF
770 GO TO 40
771 END IF
772 dif( 1 ) = dscale / dif( 1 )
773*
774* 1-norm-based estimate of Difl.
775*
776 50 CONTINUE
777 CALL dlacn2( mn2, work( mn2+1 ), work, iwork, dif( 2 ),
778 $ kase, isave )
779 IF( kase.NE.0 ) THEN
780 IF( kase.EQ.1 ) THEN
781*
782* Solve generalized Sylvester equation.
783*
784 CALL dtgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a,
785 $ lda,
786 $ work, n2, b( i, i ), ldb, b, ldb,
787 $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
788 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
789 $ ierr )
790 ELSE
791*
792* Solve the transposed variant.
793*
794 CALL dtgsyl( 'T', ijb, n2, n1, a( i, i ), lda, a,
795 $ lda,
796 $ work, n2, b( i, i ), ldb, b, ldb,
797 $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
798 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
799 $ ierr )
800 END IF
801 GO TO 50
802 END IF
803 dif( 2 ) = dscale / dif( 2 )
804*
805 END IF
806 END IF
807*
808 60 CONTINUE
809*
810* Compute generalized eigenvalues of reordered pair (A, B) and
811* normalize the generalized Schur form.
812*
813 pair = .false.
814 DO 80 k = 1, n
815 IF( pair ) THEN
816 pair = .false.
817 ELSE
818*
819 IF( k.LT.n ) THEN
820 IF( a( k+1, k ).NE.zero ) THEN
821 pair = .true.
822 END IF
823 END IF
824*
825 IF( pair ) THEN
826*
827* Compute the eigenvalue(s) at position K.
828*
829 work( 1 ) = a( k, k )
830 work( 2 ) = a( k+1, k )
831 work( 3 ) = a( k, k+1 )
832 work( 4 ) = a( k+1, k+1 )
833 work( 5 ) = b( k, k )
834 work( 6 ) = b( k+1, k )
835 work( 7 ) = b( k, k+1 )
836 work( 8 ) = b( k+1, k+1 )
837 CALL dlag2( work, 2, work( 5 ), 2, smlnum*eps,
838 $ beta( k ),
839 $ beta( k+1 ), alphar( k ), alphar( k+1 ),
840 $ alphai( k ) )
841 alphai( k+1 ) = -alphai( k )
842*
843 ELSE
844*
845 IF( sign( one, b( k, k ) ).LT.zero ) THEN
846*
847* If B(K,K) is negative, make it positive
848*
849 DO 70 i = 1, n
850 a( k, i ) = -a( k, i )
851 b( k, i ) = -b( k, i )
852 IF( wantq ) q( i, k ) = -q( i, k )
853 70 CONTINUE
854 END IF
855*
856 alphar( k ) = a( k, k )
857 alphai( k ) = zero
858 beta( k ) = b( k, k )
859*
860 END IF
861 END IF
862 80 CONTINUE
863*
864 work( 1 ) = lwmin
865 iwork( 1 ) = liwmin
866*
867 RETURN
868*
869* End of DTGSEN
870*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:134
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
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:154
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine dlassq(n, x, incx, scale, sumsq)
DLASSQ updates a sum of squares represented in scaled form.
Definition dlassq.f90:122
subroutine dtgexc(wantq, wantz, n, a, lda, b, ldb, q, ldq, z, ldz, ifst, ilst, work, lwork, info)
DTGEXC
Definition dtgexc.f:218
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:298
Here is the call graph for this function:
Here is the caller graph for this function: