Annotation of rpl/lapack/lapack/dtgex2.f, revision 1.5

1.1       bertrand    1:       SUBROUTINE DTGEX2( WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z,
                      2:      $                   LDZ, J1, N1, N2, WORK, LWORK, INFO )
                      3: *
1.5     ! bertrand    4: *  -- LAPACK auxiliary routine (version 3.2.2) --
1.1       bertrand    5: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                      6: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.5     ! bertrand    7: *     June 2010
1.1       bertrand    8: *
                      9: *     .. Scalar Arguments ..
                     10:       LOGICAL            WANTQ, WANTZ
                     11:       INTEGER            INFO, J1, LDA, LDB, LDQ, LDZ, LWORK, N, N1, N2
                     12: *     ..
                     13: *     .. Array Arguments ..
                     14:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
                     15:      $                   WORK( * ), Z( LDZ, * )
                     16: *     ..
                     17: *
                     18: *  Purpose
                     19: *  =======
                     20: *
                     21: *  DTGEX2 swaps adjacent diagonal blocks (A11, B11) and (A22, B22)
                     22: *  of size 1-by-1 or 2-by-2 in an upper (quasi) triangular matrix pair
                     23: *  (A, B) by an orthogonal equivalence transformation.
                     24: *
                     25: *  (A, B) must be in generalized real Schur canonical form (as returned
                     26: *  by DGGES), i.e. A is block upper triangular with 1-by-1 and 2-by-2
                     27: *  diagonal blocks. B is upper triangular.
                     28: *
                     29: *  Optionally, the matrices Q and Z of generalized Schur vectors are
                     30: *  updated.
                     31: *
                     32: *         Q(in) * A(in) * Z(in)' = Q(out) * A(out) * Z(out)'
                     33: *         Q(in) * B(in) * Z(in)' = Q(out) * B(out) * Z(out)'
                     34: *
                     35: *
                     36: *  Arguments
                     37: *  =========
                     38: *
                     39: *  WANTQ   (input) LOGICAL
                     40: *          .TRUE. : update the left transformation matrix Q;
                     41: *          .FALSE.: do not update Q.
                     42: *
                     43: *  WANTZ   (input) LOGICAL
                     44: *          .TRUE. : update the right transformation matrix Z;
                     45: *          .FALSE.: do not update Z.
                     46: *
                     47: *  N       (input) INTEGER
                     48: *          The order of the matrices A and B. N >= 0.
                     49: *
1.5     ! bertrand   50: *  A       (input/output) DOUBLE PRECISION array, dimensions (LDA,N)
1.1       bertrand   51: *          On entry, the matrix A in the pair (A, B).
                     52: *          On exit, the updated matrix A.
                     53: *
1.5     ! bertrand   54: *  LDA     (input) INTEGER
1.1       bertrand   55: *          The leading dimension of the array A. LDA >= max(1,N).
                     56: *
1.5     ! bertrand   57: *  B       (input/output) DOUBLE PRECISION array, dimensions (LDB,N)
1.1       bertrand   58: *          On entry, the matrix B in the pair (A, B).
                     59: *          On exit, the updated matrix B.
                     60: *
1.5     ! bertrand   61: *  LDB     (input) INTEGER
1.1       bertrand   62: *          The leading dimension of the array B. LDB >= max(1,N).
                     63: *
1.5     ! bertrand   64: *  Q       (input/output) DOUBLE PRECISION array, dimension (LDQ,N)
1.1       bertrand   65: *          On entry, if WANTQ = .TRUE., the orthogonal matrix Q.
                     66: *          On exit, the updated matrix Q.
                     67: *          Not referenced if WANTQ = .FALSE..
                     68: *
                     69: *  LDQ     (input) INTEGER
                     70: *          The leading dimension of the array Q. LDQ >= 1.
                     71: *          If WANTQ = .TRUE., LDQ >= N.
                     72: *
                     73: *  Z       (input/output) DOUBLE PRECISION array, dimension (LDZ,N)
                     74: *          On entry, if WANTZ =.TRUE., the orthogonal matrix Z.
                     75: *          On exit, the updated matrix Z.
                     76: *          Not referenced if WANTZ = .FALSE..
                     77: *
                     78: *  LDZ     (input) INTEGER
                     79: *          The leading dimension of the array Z. LDZ >= 1.
                     80: *          If WANTZ = .TRUE., LDZ >= N.
                     81: *
                     82: *  J1      (input) INTEGER
                     83: *          The index to the first block (A11, B11). 1 <= J1 <= N.
                     84: *
                     85: *  N1      (input) INTEGER
                     86: *          The order of the first block (A11, B11). N1 = 0, 1 or 2.
                     87: *
                     88: *  N2      (input) INTEGER
                     89: *          The order of the second block (A22, B22). N2 = 0, 1 or 2.
                     90: *
                     91: *  WORK    (workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)).
                     92: *
                     93: *  LWORK   (input) INTEGER
                     94: *          The dimension of the array WORK.
                     95: *          LWORK >=  MAX( 1, N*(N2+N1), (N2+N1)*(N2+N1)*2 )
                     96: *
                     97: *  INFO    (output) INTEGER
                     98: *            =0: Successful exit
                     99: *            >0: If INFO = 1, the transformed matrix (A, B) would be
                    100: *                too far from generalized Schur form; the blocks are
                    101: *                not swapped and (A, B) and (Q, Z) are unchanged.
                    102: *                The problem of swapping is too ill-conditioned.
                    103: *            <0: If INFO = -16: LWORK is too small. Appropriate value
                    104: *                for LWORK is returned in WORK(1).
                    105: *
                    106: *  Further Details
                    107: *  ===============
                    108: *
                    109: *  Based on contributions by
                    110: *     Bo Kagstrom and Peter Poromaa, Department of Computing Science,
                    111: *     Umea University, S-901 87 Umea, Sweden.
                    112: *
                    113: *  In the current code both weak and strong stability tests are
                    114: *  performed. The user can omit the strong stability test by changing
                    115: *  the internal logical parameter WANDS to .FALSE.. See ref. [2] for
                    116: *  details.
                    117: *
                    118: *  [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
                    119: *      Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
                    120: *      M.S. Moonen et al (eds), Linear Algebra for Large Scale and
                    121: *      Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
                    122: *
                    123: *  [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
                    124: *      Eigenvalues of a Regular Matrix Pair (A, B) and Condition
                    125: *      Estimation: Theory, Algorithms and Software,
                    126: *      Report UMINF - 94.04, Department of Computing Science, Umea
                    127: *      University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
                    128: *      Note 87. To appear in Numerical Algorithms, 1996.
                    129: *
                    130: *  =====================================================================
                    131: *  Replaced various illegal calls to DCOPY by calls to DLASET, or by DO
                    132: *  loops. Sven Hammarling, 1/5/02.
                    133: *
                    134: *     .. Parameters ..
                    135:       DOUBLE PRECISION   ZERO, ONE
                    136:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
1.5     ! bertrand  137:       DOUBLE PRECISION   TWENTY
        !           138:       PARAMETER          ( TWENTY = 2.0D+01 )
1.1       bertrand  139:       INTEGER            LDST
                    140:       PARAMETER          ( LDST = 4 )
                    141:       LOGICAL            WANDS
                    142:       PARAMETER          ( WANDS = .TRUE. )
                    143: *     ..
                    144: *     .. Local Scalars ..
                    145:       LOGICAL            DTRONG, WEAK
                    146:       INTEGER            I, IDUM, LINFO, M
                    147:       DOUBLE PRECISION   BQRA21, BRQA21, DDUM, DNORM, DSCALE, DSUM, EPS,
                    148:      $                   F, G, SA, SB, SCALE, SMLNUM, SS, THRESH, WS
                    149: *     ..
                    150: *     .. Local Arrays ..
                    151:       INTEGER            IWORK( LDST )
                    152:       DOUBLE PRECISION   AI( 2 ), AR( 2 ), BE( 2 ), IR( LDST, LDST ),
                    153:      $                   IRCOP( LDST, LDST ), LI( LDST, LDST ),
                    154:      $                   LICOP( LDST, LDST ), S( LDST, LDST ),
                    155:      $                   SCPY( LDST, LDST ), T( LDST, LDST ),
                    156:      $                   TAUL( LDST ), TAUR( LDST ), TCPY( LDST, LDST )
                    157: *     ..
                    158: *     .. External Functions ..
                    159:       DOUBLE PRECISION   DLAMCH
                    160:       EXTERNAL           DLAMCH
                    161: *     ..
                    162: *     .. External Subroutines ..
                    163:       EXTERNAL           DGEMM, DGEQR2, DGERQ2, DLACPY, DLAGV2, DLARTG,
                    164:      $                   DLASET, DLASSQ, DORG2R, DORGR2, DORM2R, DORMR2,
                    165:      $                   DROT, DSCAL, DTGSY2
                    166: *     ..
                    167: *     .. Intrinsic Functions ..
                    168:       INTRINSIC          ABS, MAX, SQRT
                    169: *     ..
                    170: *     .. Executable Statements ..
                    171: *
                    172:       INFO = 0
                    173: *
                    174: *     Quick return if possible
                    175: *
                    176:       IF( N.LE.1 .OR. N1.LE.0 .OR. N2.LE.0 )
                    177:      $   RETURN
                    178:       IF( N1.GT.N .OR. ( J1+N1 ).GT.N )
                    179:      $   RETURN
                    180:       M = N1 + N2
                    181:       IF( LWORK.LT.MAX( 1, N*M, M*M*2 ) ) THEN
                    182:          INFO = -16
                    183:          WORK( 1 ) = MAX( 1, N*M, M*M*2 )
                    184:          RETURN
                    185:       END IF
                    186: *
                    187:       WEAK = .FALSE.
                    188:       DTRONG = .FALSE.
                    189: *
                    190: *     Make a local copy of selected block
                    191: *
                    192:       CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, LI, LDST )
                    193:       CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, IR, LDST )
                    194:       CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, S, LDST )
                    195:       CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, T, LDST )
                    196: *
                    197: *     Compute threshold for testing acceptance of swapping.
                    198: *
                    199:       EPS = DLAMCH( 'P' )
                    200:       SMLNUM = DLAMCH( 'S' ) / EPS
                    201:       DSCALE = ZERO
                    202:       DSUM = ONE
                    203:       CALL DLACPY( 'Full', M, M, S, LDST, WORK, M )
                    204:       CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
                    205:       CALL DLACPY( 'Full', M, M, T, LDST, WORK, M )
                    206:       CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
                    207:       DNORM = DSCALE*SQRT( DSUM )
1.5     ! bertrand  208: *
        !           209: *     THRES has been changed from 
        !           210: *        THRESH = MAX( TEN*EPS*SA, SMLNUM )
        !           211: *     to
        !           212: *        THRESH = MAX( TWENTY*EPS*SA, SMLNUM )
        !           213: *     on 04/01/10.
        !           214: *     "Bug" reported by Ondra Kamenik, confirmed by Julie Langou, fixed by
        !           215: *     Jim Demmel and Guillaume Revy. See forum post 1783.
        !           216: *
        !           217:       THRESH = MAX( TWENTY*EPS*DNORM, SMLNUM )
1.1       bertrand  218: *
                    219:       IF( M.EQ.2 ) THEN
                    220: *
                    221: *        CASE 1: Swap 1-by-1 and 1-by-1 blocks.
                    222: *
                    223: *        Compute orthogonal QL and RQ that swap 1-by-1 and 1-by-1 blocks
                    224: *        using Givens rotations and perform the swap tentatively.
                    225: *
                    226:          F = S( 2, 2 )*T( 1, 1 ) - T( 2, 2 )*S( 1, 1 )
                    227:          G = S( 2, 2 )*T( 1, 2 ) - T( 2, 2 )*S( 1, 2 )
                    228:          SB = ABS( T( 2, 2 ) )
                    229:          SA = ABS( S( 2, 2 ) )
                    230:          CALL DLARTG( F, G, IR( 1, 2 ), IR( 1, 1 ), DDUM )
                    231:          IR( 2, 1 ) = -IR( 1, 2 )
                    232:          IR( 2, 2 ) = IR( 1, 1 )
                    233:          CALL DROT( 2, S( 1, 1 ), 1, S( 1, 2 ), 1, IR( 1, 1 ),
                    234:      $              IR( 2, 1 ) )
                    235:          CALL DROT( 2, T( 1, 1 ), 1, T( 1, 2 ), 1, IR( 1, 1 ),
                    236:      $              IR( 2, 1 ) )
                    237:          IF( SA.GE.SB ) THEN
                    238:             CALL DLARTG( S( 1, 1 ), S( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
                    239:      $                   DDUM )
                    240:          ELSE
                    241:             CALL DLARTG( T( 1, 1 ), T( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
                    242:      $                   DDUM )
                    243:          END IF
                    244:          CALL DROT( 2, S( 1, 1 ), LDST, S( 2, 1 ), LDST, LI( 1, 1 ),
                    245:      $              LI( 2, 1 ) )
                    246:          CALL DROT( 2, T( 1, 1 ), LDST, T( 2, 1 ), LDST, LI( 1, 1 ),
                    247:      $              LI( 2, 1 ) )
                    248:          LI( 2, 2 ) = LI( 1, 1 )
                    249:          LI( 1, 2 ) = -LI( 2, 1 )
                    250: *
                    251: *        Weak stability test:
                    252: *           |S21| + |T21| <= O(EPS * F-norm((S, T)))
                    253: *
                    254:          WS = ABS( S( 2, 1 ) ) + ABS( T( 2, 1 ) )
                    255:          WEAK = WS.LE.THRESH
                    256:          IF( .NOT.WEAK )
                    257:      $      GO TO 70
                    258: *
                    259:          IF( WANDS ) THEN
                    260: *
                    261: *           Strong stability test:
                    262: *             F-norm((A-QL'*S*QR, B-QL'*T*QR)) <= O(EPS*F-norm((A,B)))
                    263: *
                    264:             CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
                    265:      $                   M )
                    266:             CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
                    267:      $                  WORK, M )
                    268:             CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
                    269:      $                  WORK( M*M+1 ), M )
                    270:             DSCALE = ZERO
                    271:             DSUM = ONE
                    272:             CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
                    273: *
                    274:             CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
                    275:      $                   M )
                    276:             CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
                    277:      $                  WORK, M )
                    278:             CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
                    279:      $                  WORK( M*M+1 ), M )
                    280:             CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
                    281:             SS = DSCALE*SQRT( DSUM )
                    282:             DTRONG = SS.LE.THRESH
                    283:             IF( .NOT.DTRONG )
                    284:      $         GO TO 70
                    285:          END IF
                    286: *
                    287: *        Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
                    288: *               (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
                    289: *
                    290:          CALL DROT( J1+1, A( 1, J1 ), 1, A( 1, J1+1 ), 1, IR( 1, 1 ),
                    291:      $              IR( 2, 1 ) )
                    292:          CALL DROT( J1+1, B( 1, J1 ), 1, B( 1, J1+1 ), 1, IR( 1, 1 ),
                    293:      $              IR( 2, 1 ) )
                    294:          CALL DROT( N-J1+1, A( J1, J1 ), LDA, A( J1+1, J1 ), LDA,
                    295:      $              LI( 1, 1 ), LI( 2, 1 ) )
                    296:          CALL DROT( N-J1+1, B( J1, J1 ), LDB, B( J1+1, J1 ), LDB,
                    297:      $              LI( 1, 1 ), LI( 2, 1 ) )
                    298: *
                    299: *        Set  N1-by-N2 (2,1) - blocks to ZERO.
                    300: *
                    301:          A( J1+1, J1 ) = ZERO
                    302:          B( J1+1, J1 ) = ZERO
                    303: *
                    304: *        Accumulate transformations into Q and Z if requested.
                    305: *
                    306:          IF( WANTZ )
                    307:      $      CALL DROT( N, Z( 1, J1 ), 1, Z( 1, J1+1 ), 1, IR( 1, 1 ),
                    308:      $                 IR( 2, 1 ) )
                    309:          IF( WANTQ )
                    310:      $      CALL DROT( N, Q( 1, J1 ), 1, Q( 1, J1+1 ), 1, LI( 1, 1 ),
                    311:      $                 LI( 2, 1 ) )
                    312: *
                    313: *        Exit with INFO = 0 if swap was successfully performed.
                    314: *
                    315:          RETURN
                    316: *
                    317:       ELSE
                    318: *
                    319: *        CASE 2: Swap 1-by-1 and 2-by-2 blocks, or 2-by-2
                    320: *                and 2-by-2 blocks.
                    321: *
                    322: *        Solve the generalized Sylvester equation
                    323: *                 S11 * R - L * S22 = SCALE * S12
                    324: *                 T11 * R - L * T22 = SCALE * T12
                    325: *        for R and L. Solutions in LI and IR.
                    326: *
                    327:          CALL DLACPY( 'Full', N1, N2, T( 1, N1+1 ), LDST, LI, LDST )
                    328:          CALL DLACPY( 'Full', N1, N2, S( 1, N1+1 ), LDST,
                    329:      $                IR( N2+1, N1+1 ), LDST )
                    330:          CALL DTGSY2( 'N', 0, N1, N2, S, LDST, S( N1+1, N1+1 ), LDST,
                    331:      $                IR( N2+1, N1+1 ), LDST, T, LDST, T( N1+1, N1+1 ),
                    332:      $                LDST, LI, LDST, SCALE, DSUM, DSCALE, IWORK, IDUM,
                    333:      $                LINFO )
                    334: *
                    335: *        Compute orthogonal matrix QL:
                    336: *
                    337: *                    QL' * LI = [ TL ]
                    338: *                               [ 0  ]
                    339: *        where
                    340: *                    LI =  [      -L              ]
                    341: *                          [ SCALE * identity(N2) ]
                    342: *
                    343:          DO 10 I = 1, N2
                    344:             CALL DSCAL( N1, -ONE, LI( 1, I ), 1 )
                    345:             LI( N1+I, I ) = SCALE
                    346:    10    CONTINUE
                    347:          CALL DGEQR2( M, N2, LI, LDST, TAUL, WORK, LINFO )
                    348:          IF( LINFO.NE.0 )
                    349:      $      GO TO 70
                    350:          CALL DORG2R( M, M, N2, LI, LDST, TAUL, WORK, LINFO )
                    351:          IF( LINFO.NE.0 )
                    352:      $      GO TO 70
                    353: *
                    354: *        Compute orthogonal matrix RQ:
                    355: *
                    356: *                    IR * RQ' =   [ 0  TR],
                    357: *
                    358: *         where IR = [ SCALE * identity(N1), R ]
                    359: *
                    360:          DO 20 I = 1, N1
                    361:             IR( N2+I, I ) = SCALE
                    362:    20    CONTINUE
                    363:          CALL DGERQ2( N1, M, IR( N2+1, 1 ), LDST, TAUR, WORK, LINFO )
                    364:          IF( LINFO.NE.0 )
                    365:      $      GO TO 70
                    366:          CALL DORGR2( M, M, N1, IR, LDST, TAUR, WORK, LINFO )
                    367:          IF( LINFO.NE.0 )
                    368:      $      GO TO 70
                    369: *
                    370: *        Perform the swapping tentatively:
                    371: *
                    372:          CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
                    373:      $               WORK, M )
                    374:          CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, S,
                    375:      $               LDST )
                    376:          CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
                    377:      $               WORK, M )
                    378:          CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, T,
                    379:      $               LDST )
                    380:          CALL DLACPY( 'F', M, M, S, LDST, SCPY, LDST )
                    381:          CALL DLACPY( 'F', M, M, T, LDST, TCPY, LDST )
                    382:          CALL DLACPY( 'F', M, M, IR, LDST, IRCOP, LDST )
                    383:          CALL DLACPY( 'F', M, M, LI, LDST, LICOP, LDST )
                    384: *
                    385: *        Triangularize the B-part by an RQ factorization.
                    386: *        Apply transformation (from left) to A-part, giving S.
                    387: *
                    388:          CALL DGERQ2( M, M, T, LDST, TAUR, WORK, LINFO )
                    389:          IF( LINFO.NE.0 )
                    390:      $      GO TO 70
                    391:          CALL DORMR2( 'R', 'T', M, M, M, T, LDST, TAUR, S, LDST, WORK,
                    392:      $                LINFO )
                    393:          IF( LINFO.NE.0 )
                    394:      $      GO TO 70
                    395:          CALL DORMR2( 'L', 'N', M, M, M, T, LDST, TAUR, IR, LDST, WORK,
                    396:      $                LINFO )
                    397:          IF( LINFO.NE.0 )
                    398:      $      GO TO 70
                    399: *
                    400: *        Compute F-norm(S21) in BRQA21. (T21 is 0.)
                    401: *
                    402:          DSCALE = ZERO
                    403:          DSUM = ONE
                    404:          DO 30 I = 1, N2
                    405:             CALL DLASSQ( N1, S( N2+1, I ), 1, DSCALE, DSUM )
                    406:    30    CONTINUE
                    407:          BRQA21 = DSCALE*SQRT( DSUM )
                    408: *
                    409: *        Triangularize the B-part by a QR factorization.
                    410: *        Apply transformation (from right) to A-part, giving S.
                    411: *
                    412:          CALL DGEQR2( M, M, TCPY, LDST, TAUL, WORK, LINFO )
                    413:          IF( LINFO.NE.0 )
                    414:      $      GO TO 70
                    415:          CALL DORM2R( 'L', 'T', M, M, M, TCPY, LDST, TAUL, SCPY, LDST,
                    416:      $                WORK, INFO )
                    417:          CALL DORM2R( 'R', 'N', M, M, M, TCPY, LDST, TAUL, LICOP, LDST,
                    418:      $                WORK, INFO )
                    419:          IF( LINFO.NE.0 )
                    420:      $      GO TO 70
                    421: *
                    422: *        Compute F-norm(S21) in BQRA21. (T21 is 0.)
                    423: *
                    424:          DSCALE = ZERO
                    425:          DSUM = ONE
                    426:          DO 40 I = 1, N2
                    427:             CALL DLASSQ( N1, SCPY( N2+1, I ), 1, DSCALE, DSUM )
                    428:    40    CONTINUE
                    429:          BQRA21 = DSCALE*SQRT( DSUM )
                    430: *
                    431: *        Decide which method to use.
                    432: *          Weak stability test:
                    433: *             F-norm(S21) <= O(EPS * F-norm((S, T)))
                    434: *
                    435:          IF( BQRA21.LE.BRQA21 .AND. BQRA21.LE.THRESH ) THEN
                    436:             CALL DLACPY( 'F', M, M, SCPY, LDST, S, LDST )
                    437:             CALL DLACPY( 'F', M, M, TCPY, LDST, T, LDST )
                    438:             CALL DLACPY( 'F', M, M, IRCOP, LDST, IR, LDST )
                    439:             CALL DLACPY( 'F', M, M, LICOP, LDST, LI, LDST )
                    440:          ELSE IF( BRQA21.GE.THRESH ) THEN
                    441:             GO TO 70
                    442:          END IF
                    443: *
                    444: *        Set lower triangle of B-part to zero
                    445: *
                    446:          CALL DLASET( 'Lower', M-1, M-1, ZERO, ZERO, T(2,1), LDST )
                    447: *
                    448:          IF( WANDS ) THEN
                    449: *
                    450: *           Strong stability test:
                    451: *              F-norm((A-QL*S*QR', B-QL*T*QR')) <= O(EPS*F-norm((A,B)))
                    452: *
                    453:             CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
                    454:      $                   M )
                    455:             CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
                    456:      $                  WORK, M )
                    457:             CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
                    458:      $                  WORK( M*M+1 ), M )
                    459:             DSCALE = ZERO
                    460:             DSUM = ONE
                    461:             CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
                    462: *
                    463:             CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
                    464:      $                   M )
                    465:             CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
                    466:      $                  WORK, M )
                    467:             CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
                    468:      $                  WORK( M*M+1 ), M )
                    469:             CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
                    470:             SS = DSCALE*SQRT( DSUM )
                    471:             DTRONG = ( SS.LE.THRESH )
                    472:             IF( .NOT.DTRONG )
                    473:      $         GO TO 70
                    474: *
                    475:          END IF
                    476: *
                    477: *        If the swap is accepted ("weakly" and "strongly"), apply the
                    478: *        transformations and set N1-by-N2 (2,1)-block to zero.
                    479: *
                    480:          CALL DLASET( 'Full', N1, N2, ZERO, ZERO, S(N2+1,1), LDST )
                    481: *
                    482: *        copy back M-by-M diagonal block starting at index J1 of (A, B)
                    483: *
                    484:          CALL DLACPY( 'F', M, M, S, LDST, A( J1, J1 ), LDA )
                    485:          CALL DLACPY( 'F', M, M, T, LDST, B( J1, J1 ), LDB )
                    486:          CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, T, LDST )
                    487: *
                    488: *        Standardize existing 2-by-2 blocks.
                    489: *
                    490:          DO 50 I = 1, M*M
                    491:             WORK(I) = ZERO
                    492:    50    CONTINUE
                    493:          WORK( 1 ) = ONE
                    494:          T( 1, 1 ) = ONE
                    495:          IDUM = LWORK - M*M - 2
                    496:          IF( N2.GT.1 ) THEN
                    497:             CALL DLAGV2( A( J1, J1 ), LDA, B( J1, J1 ), LDB, AR, AI, BE,
                    498:      $                   WORK( 1 ), WORK( 2 ), T( 1, 1 ), T( 2, 1 ) )
                    499:             WORK( M+1 ) = -WORK( 2 )
                    500:             WORK( M+2 ) = WORK( 1 )
                    501:             T( N2, N2 ) = T( 1, 1 )
                    502:             T( 1, 2 ) = -T( 2, 1 )
                    503:          END IF
                    504:          WORK( M*M ) = ONE
                    505:          T( M, M ) = ONE
                    506: *
                    507:          IF( N1.GT.1 ) THEN
                    508:             CALL DLAGV2( A( J1+N2, J1+N2 ), LDA, B( J1+N2, J1+N2 ), LDB,
                    509:      $                   TAUR, TAUL, WORK( M*M+1 ), WORK( N2*M+N2+1 ),
                    510:      $                   WORK( N2*M+N2+2 ), T( N2+1, N2+1 ),
                    511:      $                   T( M, M-1 ) )
                    512:             WORK( M*M ) = WORK( N2*M+N2+1 )
                    513:             WORK( M*M-1 ) = -WORK( N2*M+N2+2 )
                    514:             T( M, M ) = T( N2+1, N2+1 )
                    515:             T( M-1, M ) = -T( M, M-1 )
                    516:          END IF
                    517:          CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, A( J1, J1+N2 ),
                    518:      $               LDA, ZERO, WORK( M*M+1 ), N2 )
                    519:          CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, A( J1, J1+N2 ),
                    520:      $                LDA )
                    521:          CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, B( J1, J1+N2 ),
                    522:      $               LDB, ZERO, WORK( M*M+1 ), N2 )
                    523:          CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, B( J1, J1+N2 ),
                    524:      $                LDB )
                    525:          CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, WORK, M, ZERO,
                    526:      $               WORK( M*M+1 ), M )
                    527:          CALL DLACPY( 'Full', M, M, WORK( M*M+1 ), M, LI, LDST )
                    528:          CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, A( J1, J1+N2 ), LDA,
                    529:      $               T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
                    530:          CALL DLACPY( 'Full', N2, N1, WORK, N2, A( J1, J1+N2 ), LDA )
                    531:          CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, B( J1, J1+N2 ), LDB,
                    532:      $               T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
                    533:          CALL DLACPY( 'Full', N2, N1, WORK, N2, B( J1, J1+N2 ), LDB )
                    534:          CALL DGEMM( 'T', 'N', M, M, M, ONE, IR, LDST, T, LDST, ZERO,
                    535:      $               WORK, M )
                    536:          CALL DLACPY( 'Full', M, M, WORK, M, IR, LDST )
                    537: *
                    538: *        Accumulate transformations into Q and Z if requested.
                    539: *
                    540:          IF( WANTQ ) THEN
                    541:             CALL DGEMM( 'N', 'N', N, M, M, ONE, Q( 1, J1 ), LDQ, LI,
                    542:      $                  LDST, ZERO, WORK, N )
                    543:             CALL DLACPY( 'Full', N, M, WORK, N, Q( 1, J1 ), LDQ )
                    544: *
                    545:          END IF
                    546: *
                    547:          IF( WANTZ ) THEN
                    548:             CALL DGEMM( 'N', 'N', N, M, M, ONE, Z( 1, J1 ), LDZ, IR,
                    549:      $                  LDST, ZERO, WORK, N )
                    550:             CALL DLACPY( 'Full', N, M, WORK, N, Z( 1, J1 ), LDZ )
                    551: *
                    552:          END IF
                    553: *
                    554: *        Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
                    555: *                (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
                    556: *
                    557:          I = J1 + M
                    558:          IF( I.LE.N ) THEN
                    559:             CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
                    560:      $                  A( J1, I ), LDA, ZERO, WORK, M )
                    561:             CALL DLACPY( 'Full', M, N-I+1, WORK, M, A( J1, I ), LDA )
                    562:             CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
                    563:      $                  B( J1, I ), LDA, ZERO, WORK, M )
                    564:             CALL DLACPY( 'Full', M, N-I+1, WORK, M, B( J1, I ), LDB )
                    565:          END IF
                    566:          I = J1 - 1
                    567:          IF( I.GT.0 ) THEN
                    568:             CALL DGEMM( 'N', 'N', I, M, M, ONE, A( 1, J1 ), LDA, IR,
                    569:      $                  LDST, ZERO, WORK, I )
                    570:             CALL DLACPY( 'Full', I, M, WORK, I, A( 1, J1 ), LDA )
                    571:             CALL DGEMM( 'N', 'N', I, M, M, ONE, B( 1, J1 ), LDB, IR,
                    572:      $                  LDST, ZERO, WORK, I )
                    573:             CALL DLACPY( 'Full', I, M, WORK, I, B( 1, J1 ), LDB )
                    574:          END IF
                    575: *
                    576: *        Exit with INFO = 0 if swap was successfully performed.
                    577: *
                    578:          RETURN
                    579: *
                    580:       END IF
                    581: *
                    582: *     Exit with INFO = 1 if swap was rejected.
                    583: *
                    584:    70 CONTINUE
                    585: *
                    586:       INFO = 1
                    587:       RETURN
                    588: *
                    589: *     End of DTGEX2
                    590: *
                    591:       END

CVSweb interface <joel.bertrand@systella.fr>