LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlaed3 ( integer  K,
integer  N,
integer  N1,
double precision, dimension( * )  D,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
double precision  RHO,
double precision, dimension( * )  DLAMDA,
double precision, dimension( * )  Q2,
integer, dimension( * )  INDX,
integer, dimension( * )  CTOT,
double precision, dimension( * )  W,
double precision, dimension( * )  S,
integer  INFO 
)

DLAED3 used by sstedc. Finds the roots of the secular equation and updates the eigenvectors. Used when the original matrix is tridiagonal.

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

Purpose:
 DLAED3 finds the roots of the secular equation, as defined by the
 values in D, W, and RHO, between 1 and K.  It makes the
 appropriate calls to DLAED4 and then updates the eigenvectors by
 multiplying the matrix of eigenvectors of the pair of eigensystems
 being combined by the matrix of eigenvectors of the K-by-K system
 which is solved here.

 This code makes very mild assumptions about floating point
 arithmetic. It will work on machines with a guard digit in
 add/subtract, or on those binary machines without guard digits
 which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2.
 It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[in]K
          K is INTEGER
          The number of terms in the rational function to be solved by
          DLAED4.  K >= 0.
[in]N
          N is INTEGER
          The number of rows and columns in the Q matrix.
          N >= K (deflation may result in N>K).
[in]N1
          N1 is INTEGER
          The location of the last eigenvalue in the leading submatrix.
          min(1,N) <= N1 <= N/2.
[out]D
          D is DOUBLE PRECISION array, dimension (N)
          D(I) contains the updated eigenvalues for
          1 <= I <= K.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDQ,N)
          Initially the first K columns are used as workspace.
          On output the columns 1 to K contain
          the updated eigenvectors.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  LDQ >= max(1,N).
[in]RHO
          RHO is DOUBLE PRECISION
          The value of the parameter in the rank one update equation.
          RHO >= 0 required.
[in,out]DLAMDA
          DLAMDA is DOUBLE PRECISION array, dimension (K)
          The first K elements of this array contain the old roots
          of the deflated updating problem.  These are the poles
          of the secular equation. May be changed on output by
          having lowest order bit set to zero on Cray X-MP, Cray Y-MP,
          Cray-2, or Cray C-90, as described above.
[in]Q2
          Q2 is DOUBLE PRECISION array, dimension (LDQ2, N)
          The first K columns of this matrix contain the non-deflated
          eigenvectors for the split problem.
[in]INDX
          INDX is INTEGER array, dimension (N)
          The permutation used to arrange the columns of the deflated
          Q matrix into three groups (see DLAED2).
          The rows of the eigenvectors found by DLAED4 must be likewise
          permuted before the matrix multiply can take place.
[in]CTOT
          CTOT is INTEGER array, dimension (4)
          A count of the total number of the various types of columns
          in Q, as described in INDX.  The fourth column type is any
          column which has been deflated.
[in,out]W
          W is DOUBLE PRECISION array, dimension (K)
          The first K elements of this array contain the components
          of the deflation-adjusted updating vector. Destroyed on
          output.
[out]S
          S is DOUBLE PRECISION array, dimension (N1 + 1)*K
          Will contain the eigenvectors of the repaired matrix which
          will be multiplied by the previously accumulated eigenvectors
          to update the system.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = 1, an eigenvalue did not converge
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA
Modified by Francoise Tisseur, University of Tennessee

Definition at line 187 of file dlaed3.f.

187 *
188 * -- LAPACK computational routine (version 3.4.2) --
189 * -- LAPACK is a software package provided by Univ. of Tennessee, --
190 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191 * September 2012
192 *
193 * .. Scalar Arguments ..
194  INTEGER info, k, ldq, n, n1
195  DOUBLE PRECISION rho
196 * ..
197 * .. Array Arguments ..
198  INTEGER ctot( * ), indx( * )
199  DOUBLE PRECISION d( * ), dlamda( * ), q( ldq, * ), q2( * ),
200  $ s( * ), w( * )
201 * ..
202 *
203 * =====================================================================
204 *
205 * .. Parameters ..
206  DOUBLE PRECISION one, zero
207  parameter ( one = 1.0d0, zero = 0.0d0 )
208 * ..
209 * .. Local Scalars ..
210  INTEGER i, ii, iq2, j, n12, n2, n23
211  DOUBLE PRECISION temp
212 * ..
213 * .. External Functions ..
214  DOUBLE PRECISION dlamc3, dnrm2
215  EXTERNAL dlamc3, dnrm2
216 * ..
217 * .. External Subroutines ..
218  EXTERNAL dcopy, dgemm, dlacpy, dlaed4, dlaset, xerbla
219 * ..
220 * .. Intrinsic Functions ..
221  INTRINSIC max, sign, sqrt
222 * ..
223 * .. Executable Statements ..
224 *
225 * Test the input parameters.
226 *
227  info = 0
228 *
229  IF( k.LT.0 ) THEN
230  info = -1
231  ELSE IF( n.LT.k ) THEN
232  info = -2
233  ELSE IF( ldq.LT.max( 1, n ) ) THEN
234  info = -6
235  END IF
236  IF( info.NE.0 ) THEN
237  CALL xerbla( 'DLAED3', -info )
238  RETURN
239  END IF
240 *
241 * Quick return if possible
242 *
243  IF( k.EQ.0 )
244  $ RETURN
245 *
246 * Modify values DLAMDA(i) to make sure all DLAMDA(i)-DLAMDA(j) can
247 * be computed with high relative accuracy (barring over/underflow).
248 * This is a problem on machines without a guard digit in
249 * add/subtract (Cray XMP, Cray YMP, Cray C 90 and Cray 2).
250 * The following code replaces DLAMDA(I) by 2*DLAMDA(I)-DLAMDA(I),
251 * which on any of these machines zeros out the bottommost
252 * bit of DLAMDA(I) if it is 1; this makes the subsequent
253 * subtractions DLAMDA(I)-DLAMDA(J) unproblematic when cancellation
254 * occurs. On binary machines with a guard digit (almost all
255 * machines) it does not change DLAMDA(I) at all. On hexadecimal
256 * and decimal machines with a guard digit, it slightly
257 * changes the bottommost bits of DLAMDA(I). It does not account
258 * for hexadecimal or decimal machines without guard digits
259 * (we know of none). We use a subroutine call to compute
260 * 2*DLAMBDA(I) to prevent optimizing compilers from eliminating
261 * this code.
262 *
263  DO 10 i = 1, k
264  dlamda( i ) = dlamc3( dlamda( i ), dlamda( i ) ) - dlamda( i )
265  10 CONTINUE
266 *
267  DO 20 j = 1, k
268  CALL dlaed4( k, j, dlamda, w, q( 1, j ), rho, d( j ), info )
269 *
270 * If the zero finder fails, the computation is terminated.
271 *
272  IF( info.NE.0 )
273  $ GO TO 120
274  20 CONTINUE
275 *
276  IF( k.EQ.1 )
277  $ GO TO 110
278  IF( k.EQ.2 ) THEN
279  DO 30 j = 1, k
280  w( 1 ) = q( 1, j )
281  w( 2 ) = q( 2, j )
282  ii = indx( 1 )
283  q( 1, j ) = w( ii )
284  ii = indx( 2 )
285  q( 2, j ) = w( ii )
286  30 CONTINUE
287  GO TO 110
288  END IF
289 *
290 * Compute updated W.
291 *
292  CALL dcopy( k, w, 1, s, 1 )
293 *
294 * Initialize W(I) = Q(I,I)
295 *
296  CALL dcopy( k, q, ldq+1, w, 1 )
297  DO 60 j = 1, k
298  DO 40 i = 1, j - 1
299  w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
300  40 CONTINUE
301  DO 50 i = j + 1, k
302  w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
303  50 CONTINUE
304  60 CONTINUE
305  DO 70 i = 1, k
306  w( i ) = sign( sqrt( -w( i ) ), s( i ) )
307  70 CONTINUE
308 *
309 * Compute eigenvectors of the modified rank-1 modification.
310 *
311  DO 100 j = 1, k
312  DO 80 i = 1, k
313  s( i ) = w( i ) / q( i, j )
314  80 CONTINUE
315  temp = dnrm2( k, s, 1 )
316  DO 90 i = 1, k
317  ii = indx( i )
318  q( i, j ) = s( ii ) / temp
319  90 CONTINUE
320  100 CONTINUE
321 *
322 * Compute the updated eigenvectors.
323 *
324  110 CONTINUE
325 *
326  n2 = n - n1
327  n12 = ctot( 1 ) + ctot( 2 )
328  n23 = ctot( 2 ) + ctot( 3 )
329 *
330  CALL dlacpy( 'A', n23, k, q( ctot( 1 )+1, 1 ), ldq, s, n23 )
331  iq2 = n1*n12 + 1
332  IF( n23.NE.0 ) THEN
333  CALL dgemm( 'N', 'N', n2, k, n23, one, q2( iq2 ), n2, s, n23,
334  $ zero, q( n1+1, 1 ), ldq )
335  ELSE
336  CALL dlaset( 'A', n2, k, zero, zero, q( n1+1, 1 ), ldq )
337  END IF
338 *
339  CALL dlacpy( 'A', n12, k, q, ldq, s, n12 )
340  IF( n12.NE.0 ) THEN
341  CALL dgemm( 'N', 'N', n1, k, n12, one, q2, n1, s, n12, zero, q,
342  $ ldq )
343  ELSE
344  CALL dlaset( 'A', n1, k, zero, zero, q( 1, 1 ), ldq )
345  END IF
346 *
347 *
348  120 CONTINUE
349  RETURN
350 *
351 * End of DLAED3
352 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
double precision function dlamc3(A, B)
DLAMC3
Definition: dlamch.f:169
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
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 dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
double precision function dnrm2(N, X, INCX)
DNRM2
Definition: dnrm2.f:56
subroutine dlaed4(N, I, D, Z, DELTA, RHO, DLAM, INFO)
DLAED4 used by sstedc. Finds a single root of the secular equation.
Definition: dlaed4.f:147

Here is the call graph for this function:

Here is the caller graph for this function: