Annotation of rpl/lapack/lapack/zlaed8.f, revision 1.4

1.1       bertrand    1:       SUBROUTINE ZLAED8( K, N, QSIZ, Q, LDQ, D, RHO, CUTPNT, Z, DLAMDA,
                      2:      $                   Q2, LDQ2, W, INDXP, INDX, INDXQ, PERM, GIVPTR,
                      3:      $                   GIVCOL, GIVNUM, INFO )
                      4: *
                      5: *  -- LAPACK routine (version 3.2) --
                      6: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                      7: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
                      8: *     November 2006
                      9: *
                     10: *     .. Scalar Arguments ..
                     11:       INTEGER            CUTPNT, GIVPTR, INFO, K, LDQ, LDQ2, N, QSIZ
                     12:       DOUBLE PRECISION   RHO
                     13: *     ..
                     14: *     .. Array Arguments ..
                     15:       INTEGER            GIVCOL( 2, * ), INDX( * ), INDXP( * ),
                     16:      $                   INDXQ( * ), PERM( * )
                     17:       DOUBLE PRECISION   D( * ), DLAMDA( * ), GIVNUM( 2, * ), W( * ),
                     18:      $                   Z( * )
                     19:       COMPLEX*16         Q( LDQ, * ), Q2( LDQ2, * )
                     20: *     ..
                     21: *
                     22: *  Purpose
                     23: *  =======
                     24: *
                     25: *  ZLAED8 merges the two sets of eigenvalues together into a single
                     26: *  sorted set.  Then it tries to deflate the size of the problem.
                     27: *  There are two ways in which deflation can occur:  when two or more
                     28: *  eigenvalues are close together or if there is a tiny element in the
                     29: *  Z vector.  For each such occurrence the order of the related secular
                     30: *  equation problem is reduced by one.
                     31: *
                     32: *  Arguments
                     33: *  =========
                     34: *
                     35: *  K      (output) INTEGER
                     36: *         Contains the number of non-deflated eigenvalues.
                     37: *         This is the order of the related secular equation.
                     38: *
                     39: *  N      (input) INTEGER
                     40: *         The dimension of the symmetric tridiagonal matrix.  N >= 0.
                     41: *
                     42: *  QSIZ   (input) INTEGER
                     43: *         The dimension of the unitary matrix used to reduce
                     44: *         the dense or band matrix to tridiagonal form.
                     45: *         QSIZ >= N if ICOMPQ = 1.
                     46: *
                     47: *  Q      (input/output) COMPLEX*16 array, dimension (LDQ,N)
                     48: *         On entry, Q contains the eigenvectors of the partially solved
                     49: *         system which has been previously updated in matrix
                     50: *         multiplies with other partially solved eigensystems.
                     51: *         On exit, Q contains the trailing (N-K) updated eigenvectors
                     52: *         (those which were deflated) in its last N-K columns.
                     53: *
                     54: *  LDQ    (input) INTEGER
                     55: *         The leading dimension of the array Q.  LDQ >= max( 1, N ).
                     56: *
                     57: *  D      (input/output) DOUBLE PRECISION array, dimension (N)
                     58: *         On entry, D contains the eigenvalues of the two submatrices to
                     59: *         be combined.  On exit, D contains the trailing (N-K) updated
                     60: *         eigenvalues (those which were deflated) sorted into increasing
                     61: *         order.
                     62: *
                     63: *  RHO    (input/output) DOUBLE PRECISION
                     64: *         Contains the off diagonal element associated with the rank-1
                     65: *         cut which originally split the two submatrices which are now
                     66: *         being recombined. RHO is modified during the computation to
                     67: *         the value required by DLAED3.
                     68: *
                     69: *  CUTPNT (input) INTEGER
                     70: *         Contains the location of the last eigenvalue in the leading
                     71: *         sub-matrix.  MIN(1,N) <= CUTPNT <= N.
                     72: *
                     73: *  Z      (input) DOUBLE PRECISION array, dimension (N)
                     74: *         On input this vector contains the updating vector (the last
                     75: *         row of the first sub-eigenvector matrix and the first row of
                     76: *         the second sub-eigenvector matrix).  The contents of Z are
                     77: *         destroyed during the updating process.
                     78: *
                     79: *  DLAMDA (output) DOUBLE PRECISION array, dimension (N)
                     80: *         Contains a copy of the first K eigenvalues which will be used
                     81: *         by DLAED3 to form the secular equation.
                     82: *
                     83: *  Q2     (output) COMPLEX*16 array, dimension (LDQ2,N)
                     84: *         If ICOMPQ = 0, Q2 is not referenced.  Otherwise,
                     85: *         Contains a copy of the first K eigenvectors which will be used
                     86: *         by DLAED7 in a matrix multiply (DGEMM) to update the new
                     87: *         eigenvectors.
                     88: *
                     89: *  LDQ2   (input) INTEGER
                     90: *         The leading dimension of the array Q2.  LDQ2 >= max( 1, N ).
                     91: *
                     92: *  W      (output) DOUBLE PRECISION array, dimension (N)
                     93: *         This will hold the first k values of the final
                     94: *         deflation-altered z-vector and will be passed to DLAED3.
                     95: *
                     96: *  INDXP  (workspace) INTEGER array, dimension (N)
                     97: *         This will contain the permutation used to place deflated
                     98: *         values of D at the end of the array. On output INDXP(1:K)
                     99: *         points to the nondeflated D-values and INDXP(K+1:N)
                    100: *         points to the deflated eigenvalues.
                    101: *
                    102: *  INDX   (workspace) INTEGER array, dimension (N)
                    103: *         This will contain the permutation used to sort the contents of
                    104: *         D into ascending order.
                    105: *
                    106: *  INDXQ  (input) INTEGER array, dimension (N)
                    107: *         This contains the permutation which separately sorts the two
                    108: *         sub-problems in D into ascending order.  Note that elements in
                    109: *         the second half of this permutation must first have CUTPNT
                    110: *         added to their values in order to be accurate.
                    111: *
                    112: *  PERM   (output) INTEGER array, dimension (N)
                    113: *         Contains the permutations (from deflation and sorting) to be
                    114: *         applied to each eigenblock.
                    115: *
                    116: *  GIVPTR (output) INTEGER
                    117: *         Contains the number of Givens rotations which took place in
                    118: *         this subproblem.
                    119: *
                    120: *  GIVCOL (output) INTEGER array, dimension (2, N)
                    121: *         Each pair of numbers indicates a pair of columns to take place
                    122: *         in a Givens rotation.
                    123: *
                    124: *  GIVNUM (output) DOUBLE PRECISION array, dimension (2, N)
                    125: *         Each number indicates the S value to be used in the
                    126: *         corresponding Givens rotation.
                    127: *
                    128: *  INFO   (output) INTEGER
                    129: *          = 0:  successful exit.
                    130: *          < 0:  if INFO = -i, the i-th argument had an illegal value.
                    131: *
                    132: *  =====================================================================
                    133: *
                    134: *     .. Parameters ..
                    135:       DOUBLE PRECISION   MONE, ZERO, ONE, TWO, EIGHT
                    136:       PARAMETER          ( MONE = -1.0D0, ZERO = 0.0D0, ONE = 1.0D0,
                    137:      $                   TWO = 2.0D0, EIGHT = 8.0D0 )
                    138: *     ..
                    139: *     .. Local Scalars ..
                    140:       INTEGER            I, IMAX, J, JLAM, JMAX, JP, K2, N1, N1P1, N2
                    141:       DOUBLE PRECISION   C, EPS, S, T, TAU, TOL
                    142: *     ..
                    143: *     .. External Functions ..
                    144:       INTEGER            IDAMAX
                    145:       DOUBLE PRECISION   DLAMCH, DLAPY2
                    146:       EXTERNAL           IDAMAX, DLAMCH, DLAPY2
                    147: *     ..
                    148: *     .. External Subroutines ..
                    149:       EXTERNAL           DCOPY, DLAMRG, DSCAL, XERBLA, ZCOPY, ZDROT,
                    150:      $                   ZLACPY
                    151: *     ..
                    152: *     .. Intrinsic Functions ..
                    153:       INTRINSIC          ABS, MAX, MIN, SQRT
                    154: *     ..
                    155: *     .. Executable Statements ..
                    156: *
                    157: *     Test the input parameters.
                    158: *
                    159:       INFO = 0
                    160: *
                    161:       IF( N.LT.0 ) THEN
                    162:          INFO = -2
                    163:       ELSE IF( QSIZ.LT.N ) THEN
                    164:          INFO = -3
                    165:       ELSE IF( LDQ.LT.MAX( 1, N ) ) THEN
                    166:          INFO = -5
                    167:       ELSE IF( CUTPNT.LT.MIN( 1, N ) .OR. CUTPNT.GT.N ) THEN
                    168:          INFO = -8
                    169:       ELSE IF( LDQ2.LT.MAX( 1, N ) ) THEN
                    170:          INFO = -12
                    171:       END IF
                    172:       IF( INFO.NE.0 ) THEN
                    173:          CALL XERBLA( 'ZLAED8', -INFO )
                    174:          RETURN
                    175:       END IF
                    176: *
                    177: *     Quick return if possible
                    178: *
                    179:       IF( N.EQ.0 )
                    180:      $   RETURN
                    181: *
                    182:       N1 = CUTPNT
                    183:       N2 = N - N1
                    184:       N1P1 = N1 + 1
                    185: *
                    186:       IF( RHO.LT.ZERO ) THEN
                    187:          CALL DSCAL( N2, MONE, Z( N1P1 ), 1 )
                    188:       END IF
                    189: *
                    190: *     Normalize z so that norm(z) = 1
                    191: *
                    192:       T = ONE / SQRT( TWO )
                    193:       DO 10 J = 1, N
                    194:          INDX( J ) = J
                    195:    10 CONTINUE
                    196:       CALL DSCAL( N, T, Z, 1 )
                    197:       RHO = ABS( TWO*RHO )
                    198: *
                    199: *     Sort the eigenvalues into increasing order
                    200: *
                    201:       DO 20 I = CUTPNT + 1, N
                    202:          INDXQ( I ) = INDXQ( I ) + CUTPNT
                    203:    20 CONTINUE
                    204:       DO 30 I = 1, N
                    205:          DLAMDA( I ) = D( INDXQ( I ) )
                    206:          W( I ) = Z( INDXQ( I ) )
                    207:    30 CONTINUE
                    208:       I = 1
                    209:       J = CUTPNT + 1
                    210:       CALL DLAMRG( N1, N2, DLAMDA, 1, 1, INDX )
                    211:       DO 40 I = 1, N
                    212:          D( I ) = DLAMDA( INDX( I ) )
                    213:          Z( I ) = W( INDX( I ) )
                    214:    40 CONTINUE
                    215: *
                    216: *     Calculate the allowable deflation tolerance
                    217: *
                    218:       IMAX = IDAMAX( N, Z, 1 )
                    219:       JMAX = IDAMAX( N, D, 1 )
                    220:       EPS = DLAMCH( 'Epsilon' )
                    221:       TOL = EIGHT*EPS*ABS( D( JMAX ) )
                    222: *
                    223: *     If the rank-1 modifier is small enough, no more needs to be done
                    224: *     -- except to reorganize Q so that its columns correspond with the
                    225: *     elements in D.
                    226: *
                    227:       IF( RHO*ABS( Z( IMAX ) ).LE.TOL ) THEN
                    228:          K = 0
                    229:          DO 50 J = 1, N
                    230:             PERM( J ) = INDXQ( INDX( J ) )
                    231:             CALL ZCOPY( QSIZ, Q( 1, PERM( J ) ), 1, Q2( 1, J ), 1 )
                    232:    50    CONTINUE
                    233:          CALL ZLACPY( 'A', QSIZ, N, Q2( 1, 1 ), LDQ2, Q( 1, 1 ), LDQ )
                    234:          RETURN
                    235:       END IF
                    236: *
                    237: *     If there are multiple eigenvalues then the problem deflates.  Here
                    238: *     the number of equal eigenvalues are found.  As each equal
                    239: *     eigenvalue is found, an elementary reflector is computed to rotate
                    240: *     the corresponding eigensubspace so that the corresponding
                    241: *     components of Z are zero in this new basis.
                    242: *
                    243:       K = 0
                    244:       GIVPTR = 0
                    245:       K2 = N + 1
                    246:       DO 60 J = 1, N
                    247:          IF( RHO*ABS( Z( J ) ).LE.TOL ) THEN
                    248: *
                    249: *           Deflate due to small z component.
                    250: *
                    251:             K2 = K2 - 1
                    252:             INDXP( K2 ) = J
                    253:             IF( J.EQ.N )
                    254:      $         GO TO 100
                    255:          ELSE
                    256:             JLAM = J
                    257:             GO TO 70
                    258:          END IF
                    259:    60 CONTINUE
                    260:    70 CONTINUE
                    261:       J = J + 1
                    262:       IF( J.GT.N )
                    263:      $   GO TO 90
                    264:       IF( RHO*ABS( Z( J ) ).LE.TOL ) THEN
                    265: *
                    266: *        Deflate due to small z component.
                    267: *
                    268:          K2 = K2 - 1
                    269:          INDXP( K2 ) = J
                    270:       ELSE
                    271: *
                    272: *        Check if eigenvalues are close enough to allow deflation.
                    273: *
                    274:          S = Z( JLAM )
                    275:          C = Z( J )
                    276: *
                    277: *        Find sqrt(a**2+b**2) without overflow or
                    278: *        destructive underflow.
                    279: *
                    280:          TAU = DLAPY2( C, S )
                    281:          T = D( J ) - D( JLAM )
                    282:          C = C / TAU
                    283:          S = -S / TAU
                    284:          IF( ABS( T*C*S ).LE.TOL ) THEN
                    285: *
                    286: *           Deflation is possible.
                    287: *
                    288:             Z( J ) = TAU
                    289:             Z( JLAM ) = ZERO
                    290: *
                    291: *           Record the appropriate Givens rotation
                    292: *
                    293:             GIVPTR = GIVPTR + 1
                    294:             GIVCOL( 1, GIVPTR ) = INDXQ( INDX( JLAM ) )
                    295:             GIVCOL( 2, GIVPTR ) = INDXQ( INDX( J ) )
                    296:             GIVNUM( 1, GIVPTR ) = C
                    297:             GIVNUM( 2, GIVPTR ) = S
                    298:             CALL ZDROT( QSIZ, Q( 1, INDXQ( INDX( JLAM ) ) ), 1,
                    299:      $                  Q( 1, INDXQ( INDX( J ) ) ), 1, C, S )
                    300:             T = D( JLAM )*C*C + D( J )*S*S
                    301:             D( J ) = D( JLAM )*S*S + D( J )*C*C
                    302:             D( JLAM ) = T
                    303:             K2 = K2 - 1
                    304:             I = 1
                    305:    80       CONTINUE
                    306:             IF( K2+I.LE.N ) THEN
                    307:                IF( D( JLAM ).LT.D( INDXP( K2+I ) ) ) THEN
                    308:                   INDXP( K2+I-1 ) = INDXP( K2+I )
                    309:                   INDXP( K2+I ) = JLAM
                    310:                   I = I + 1
                    311:                   GO TO 80
                    312:                ELSE
                    313:                   INDXP( K2+I-1 ) = JLAM
                    314:                END IF
                    315:             ELSE
                    316:                INDXP( K2+I-1 ) = JLAM
                    317:             END IF
                    318:             JLAM = J
                    319:          ELSE
                    320:             K = K + 1
                    321:             W( K ) = Z( JLAM )
                    322:             DLAMDA( K ) = D( JLAM )
                    323:             INDXP( K ) = JLAM
                    324:             JLAM = J
                    325:          END IF
                    326:       END IF
                    327:       GO TO 70
                    328:    90 CONTINUE
                    329: *
                    330: *     Record the last eigenvalue.
                    331: *
                    332:       K = K + 1
                    333:       W( K ) = Z( JLAM )
                    334:       DLAMDA( K ) = D( JLAM )
                    335:       INDXP( K ) = JLAM
                    336: *
                    337:   100 CONTINUE
                    338: *
                    339: *     Sort the eigenvalues and corresponding eigenvectors into DLAMDA
                    340: *     and Q2 respectively.  The eigenvalues/vectors which were not
                    341: *     deflated go into the first K slots of DLAMDA and Q2 respectively,
                    342: *     while those which were deflated go into the last N - K slots.
                    343: *
                    344:       DO 110 J = 1, N
                    345:          JP = INDXP( J )
                    346:          DLAMDA( J ) = D( JP )
                    347:          PERM( J ) = INDXQ( INDX( JP ) )
                    348:          CALL ZCOPY( QSIZ, Q( 1, PERM( J ) ), 1, Q2( 1, J ), 1 )
                    349:   110 CONTINUE
                    350: *
                    351: *     The deflated eigenvalues and their corresponding vectors go back
                    352: *     into the last N - K slots of D and Q respectively.
                    353: *
                    354:       IF( K.LT.N ) THEN
                    355:          CALL DCOPY( N-K, DLAMDA( K+1 ), 1, D( K+1 ), 1 )
                    356:          CALL ZLACPY( 'A', QSIZ, N-K, Q2( 1, K+1 ), LDQ2, Q( 1, K+1 ),
                    357:      $                LDQ )
                    358:       END IF
                    359: *
                    360:       RETURN
                    361: *
                    362: *     End of ZLAED8
                    363: *
                    364:       END

CVSweb interface <joel.bertrand@systella.fr>