LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlaed8 ( integer  ICOMPQ,
integer  K,
integer  N,
integer  QSIZ,
double precision, dimension( * )  D,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
integer, dimension( * )  INDXQ,
double precision  RHO,
integer  CUTPNT,
double precision, dimension( * )  Z,
double precision, dimension( * )  DLAMDA,
double precision, dimension( ldq2, * )  Q2,
integer  LDQ2,
double precision, dimension( * )  W,
integer, dimension( * )  PERM,
integer  GIVPTR,
integer, dimension( 2, * )  GIVCOL,
double precision, dimension( 2, * )  GIVNUM,
integer, dimension( * )  INDXP,
integer, dimension( * )  INDX,
integer  INFO 
)

DLAED8 used by sstedc. Merges eigenvalues and deflates secular equation. Used when the original matrix is dense.

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

Purpose:
 DLAED8 merges the two sets of eigenvalues together into a single
 sorted set.  Then it tries to deflate the size of the problem.
 There are two ways in which deflation can occur:  when two or more
 eigenvalues are close together or if there is a tiny element in the
 Z vector.  For each such occurrence the order of the related secular
 equation problem is reduced by one.
Parameters
[in]ICOMPQ
          ICOMPQ is INTEGER
          = 0:  Compute eigenvalues only.
          = 1:  Compute eigenvectors of original dense symmetric matrix
                also.  On entry, Q contains the orthogonal matrix used
                to reduce the original matrix to tridiagonal form.
[out]K
          K is INTEGER
         The number of non-deflated eigenvalues, and the order of the
         related secular equation.
[in]N
          N is INTEGER
         The dimension of the symmetric tridiagonal matrix.  N >= 0.
[in]QSIZ
          QSIZ is INTEGER
         The dimension of the orthogonal matrix used to reduce
         the full matrix to tridiagonal form.  QSIZ >= N if ICOMPQ = 1.
[in,out]D
          D is DOUBLE PRECISION array, dimension (N)
         On entry, the eigenvalues of the two submatrices to be
         combined.  On exit, the trailing (N-K) updated eigenvalues
         (those which were deflated) sorted into increasing order.
[in,out]Q
          Q is DOUBLE PRECISION array, dimension (LDQ,N)
         If ICOMPQ = 0, Q is not referenced.  Otherwise,
         on entry, Q contains the eigenvectors of the partially solved
         system which has been previously updated in matrix
         multiplies with other partially solved eigensystems.
         On exit, Q contains the trailing (N-K) updated eigenvectors
         (those which were deflated) in its last N-K columns.
[in]LDQ
          LDQ is INTEGER
         The leading dimension of the array Q.  LDQ >= max(1,N).
[in]INDXQ
          INDXQ is INTEGER array, dimension (N)
         The permutation which separately sorts the two sub-problems
         in D into ascending order.  Note that elements in the second
         half of this permutation must first have CUTPNT added to
         their values in order to be accurate.
[in,out]RHO
          RHO is DOUBLE PRECISION
         On entry, the off-diagonal element associated with the rank-1
         cut which originally split the two submatrices which are now
         being recombined.
         On exit, RHO has been modified to the value required by
         DLAED3.
[in]CUTPNT
          CUTPNT is INTEGER
         The location of the last eigenvalue in the leading
         sub-matrix.  min(1,N) <= CUTPNT <= N.
[in]Z
          Z is DOUBLE PRECISION array, dimension (N)
         On entry, Z contains the updating vector (the last row of
         the first sub-eigenvector matrix and the first row of the
         second sub-eigenvector matrix).
         On exit, the contents of Z are destroyed by the updating
         process.
[out]DLAMDA
          DLAMDA is DOUBLE PRECISION array, dimension (N)
         A copy of the first K eigenvalues which will be used by
         DLAED3 to form the secular equation.
[out]Q2
          Q2 is DOUBLE PRECISION array, dimension (LDQ2,N)
         If ICOMPQ = 0, Q2 is not referenced.  Otherwise,
         a copy of the first K eigenvectors which will be used by
         DLAED7 in a matrix multiply (DGEMM) to update the new
         eigenvectors.
[in]LDQ2
          LDQ2 is INTEGER
         The leading dimension of the array Q2.  LDQ2 >= max(1,N).
[out]W
          W is DOUBLE PRECISION array, dimension (N)
         The first k values of the final deflation-altered z-vector and
         will be passed to DLAED3.
[out]PERM
          PERM is INTEGER array, dimension (N)
         The permutations (from deflation and sorting) to be applied
         to each eigenblock.
[out]GIVPTR
          GIVPTR is INTEGER
         The number of Givens rotations which took place in this
         subproblem.
[out]GIVCOL
          GIVCOL is INTEGER array, dimension (2, N)
         Each pair of numbers indicates a pair of columns to take place
         in a Givens rotation.
[out]GIVNUM
          GIVNUM is DOUBLE PRECISION array, dimension (2, N)
         Each number indicates the S value to be used in the
         corresponding Givens rotation.
[out]INDXP
          INDXP is INTEGER array, dimension (N)
         The permutation used to place deflated values of D at the end
         of the array.  INDXP(1:K) points to the nondeflated D-values
         and INDXP(K+1:N) points to the deflated eigenvalues.
[out]INDX
          INDX is INTEGER array, dimension (N)
         The permutation used to sort the contents of D into ascending
         order.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
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

Definition at line 245 of file dlaed8.f.

245 *
246 * -- LAPACK computational routine (version 3.4.2) --
247 * -- LAPACK is a software package provided by Univ. of Tennessee, --
248 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
249 * September 2012
250 *
251 * .. Scalar Arguments ..
252  INTEGER cutpnt, givptr, icompq, info, k, ldq, ldq2, n,
253  $ qsiz
254  DOUBLE PRECISION rho
255 * ..
256 * .. Array Arguments ..
257  INTEGER givcol( 2, * ), indx( * ), indxp( * ),
258  $ indxq( * ), perm( * )
259  DOUBLE PRECISION d( * ), dlamda( * ), givnum( 2, * ),
260  $ q( ldq, * ), q2( ldq2, * ), w( * ), z( * )
261 * ..
262 *
263 * =====================================================================
264 *
265 * .. Parameters ..
266  DOUBLE PRECISION mone, zero, one, two, eight
267  parameter ( mone = -1.0d0, zero = 0.0d0, one = 1.0d0,
268  $ two = 2.0d0, eight = 8.0d0 )
269 * ..
270 * .. Local Scalars ..
271 *
272  INTEGER i, imax, j, jlam, jmax, jp, k2, n1, n1p1, n2
273  DOUBLE PRECISION c, eps, s, t, tau, tol
274 * ..
275 * .. External Functions ..
276  INTEGER idamax
277  DOUBLE PRECISION dlamch, dlapy2
278  EXTERNAL idamax, dlamch, dlapy2
279 * ..
280 * .. External Subroutines ..
281  EXTERNAL dcopy, dlacpy, dlamrg, drot, dscal, xerbla
282 * ..
283 * .. Intrinsic Functions ..
284  INTRINSIC abs, max, min, sqrt
285 * ..
286 * .. Executable Statements ..
287 *
288 * Test the input parameters.
289 *
290  info = 0
291 *
292  IF( icompq.LT.0 .OR. icompq.GT.1 ) THEN
293  info = -1
294  ELSE IF( n.LT.0 ) THEN
295  info = -3
296  ELSE IF( icompq.EQ.1 .AND. qsiz.LT.n ) THEN
297  info = -4
298  ELSE IF( ldq.LT.max( 1, n ) ) THEN
299  info = -7
300  ELSE IF( cutpnt.LT.min( 1, n ) .OR. cutpnt.GT.n ) THEN
301  info = -10
302  ELSE IF( ldq2.LT.max( 1, n ) ) THEN
303  info = -14
304  END IF
305  IF( info.NE.0 ) THEN
306  CALL xerbla( 'DLAED8', -info )
307  RETURN
308  END IF
309 *
310 * Need to initialize GIVPTR to O here in case of quick exit
311 * to prevent an unspecified code behavior (usually sigfault)
312 * when IWORK array on entry to *stedc is not zeroed
313 * (or at least some IWORK entries which used in *laed7 for GIVPTR).
314 *
315  givptr = 0
316 *
317 * Quick return if possible
318 *
319  IF( n.EQ.0 )
320  $ RETURN
321 *
322  n1 = cutpnt
323  n2 = n - n1
324  n1p1 = n1 + 1
325 *
326  IF( rho.LT.zero ) THEN
327  CALL dscal( n2, mone, z( n1p1 ), 1 )
328  END IF
329 *
330 * Normalize z so that norm(z) = 1
331 *
332  t = one / sqrt( two )
333  DO 10 j = 1, n
334  indx( j ) = j
335  10 CONTINUE
336  CALL dscal( n, t, z, 1 )
337  rho = abs( two*rho )
338 *
339 * Sort the eigenvalues into increasing order
340 *
341  DO 20 i = cutpnt + 1, n
342  indxq( i ) = indxq( i ) + cutpnt
343  20 CONTINUE
344  DO 30 i = 1, n
345  dlamda( i ) = d( indxq( i ) )
346  w( i ) = z( indxq( i ) )
347  30 CONTINUE
348  i = 1
349  j = cutpnt + 1
350  CALL dlamrg( n1, n2, dlamda, 1, 1, indx )
351  DO 40 i = 1, n
352  d( i ) = dlamda( indx( i ) )
353  z( i ) = w( indx( i ) )
354  40 CONTINUE
355 *
356 * Calculate the allowable deflation tolerence
357 *
358  imax = idamax( n, z, 1 )
359  jmax = idamax( n, d, 1 )
360  eps = dlamch( 'Epsilon' )
361  tol = eight*eps*abs( d( jmax ) )
362 *
363 * If the rank-1 modifier is small enough, no more needs to be done
364 * except to reorganize Q so that its columns correspond with the
365 * elements in D.
366 *
367  IF( rho*abs( z( imax ) ).LE.tol ) THEN
368  k = 0
369  IF( icompq.EQ.0 ) THEN
370  DO 50 j = 1, n
371  perm( j ) = indxq( indx( j ) )
372  50 CONTINUE
373  ELSE
374  DO 60 j = 1, n
375  perm( j ) = indxq( indx( j ) )
376  CALL dcopy( qsiz, q( 1, perm( j ) ), 1, q2( 1, j ), 1 )
377  60 CONTINUE
378  CALL dlacpy( 'A', qsiz, n, q2( 1, 1 ), ldq2, q( 1, 1 ),
379  $ ldq )
380  END IF
381  RETURN
382  END IF
383 *
384 * If there are multiple eigenvalues then the problem deflates. Here
385 * the number of equal eigenvalues are found. As each equal
386 * eigenvalue is found, an elementary reflector is computed to rotate
387 * the corresponding eigensubspace so that the corresponding
388 * components of Z are zero in this new basis.
389 *
390  k = 0
391  k2 = n + 1
392  DO 70 j = 1, n
393  IF( rho*abs( z( j ) ).LE.tol ) THEN
394 *
395 * Deflate due to small z component.
396 *
397  k2 = k2 - 1
398  indxp( k2 ) = j
399  IF( j.EQ.n )
400  $ GO TO 110
401  ELSE
402  jlam = j
403  GO TO 80
404  END IF
405  70 CONTINUE
406  80 CONTINUE
407  j = j + 1
408  IF( j.GT.n )
409  $ GO TO 100
410  IF( rho*abs( z( j ) ).LE.tol ) THEN
411 *
412 * Deflate due to small z component.
413 *
414  k2 = k2 - 1
415  indxp( k2 ) = j
416  ELSE
417 *
418 * Check if eigenvalues are close enough to allow deflation.
419 *
420  s = z( jlam )
421  c = z( j )
422 *
423 * Find sqrt(a**2+b**2) without overflow or
424 * destructive underflow.
425 *
426  tau = dlapy2( c, s )
427  t = d( j ) - d( jlam )
428  c = c / tau
429  s = -s / tau
430  IF( abs( t*c*s ).LE.tol ) THEN
431 *
432 * Deflation is possible.
433 *
434  z( j ) = tau
435  z( jlam ) = zero
436 *
437 * Record the appropriate Givens rotation
438 *
439  givptr = givptr + 1
440  givcol( 1, givptr ) = indxq( indx( jlam ) )
441  givcol( 2, givptr ) = indxq( indx( j ) )
442  givnum( 1, givptr ) = c
443  givnum( 2, givptr ) = s
444  IF( icompq.EQ.1 ) THEN
445  CALL drot( qsiz, q( 1, indxq( indx( jlam ) ) ), 1,
446  $ q( 1, indxq( indx( j ) ) ), 1, c, s )
447  END IF
448  t = d( jlam )*c*c + d( j )*s*s
449  d( j ) = d( jlam )*s*s + d( j )*c*c
450  d( jlam ) = t
451  k2 = k2 - 1
452  i = 1
453  90 CONTINUE
454  IF( k2+i.LE.n ) THEN
455  IF( d( jlam ).LT.d( indxp( k2+i ) ) ) THEN
456  indxp( k2+i-1 ) = indxp( k2+i )
457  indxp( k2+i ) = jlam
458  i = i + 1
459  GO TO 90
460  ELSE
461  indxp( k2+i-1 ) = jlam
462  END IF
463  ELSE
464  indxp( k2+i-1 ) = jlam
465  END IF
466  jlam = j
467  ELSE
468  k = k + 1
469  w( k ) = z( jlam )
470  dlamda( k ) = d( jlam )
471  indxp( k ) = jlam
472  jlam = j
473  END IF
474  END IF
475  GO TO 80
476  100 CONTINUE
477 *
478 * Record the last eigenvalue.
479 *
480  k = k + 1
481  w( k ) = z( jlam )
482  dlamda( k ) = d( jlam )
483  indxp( k ) = jlam
484 *
485  110 CONTINUE
486 *
487 * Sort the eigenvalues and corresponding eigenvectors into DLAMDA
488 * and Q2 respectively. The eigenvalues/vectors which were not
489 * deflated go into the first K slots of DLAMDA and Q2 respectively,
490 * while those which were deflated go into the last N - K slots.
491 *
492  IF( icompq.EQ.0 ) THEN
493  DO 120 j = 1, n
494  jp = indxp( j )
495  dlamda( j ) = d( jp )
496  perm( j ) = indxq( indx( jp ) )
497  120 CONTINUE
498  ELSE
499  DO 130 j = 1, n
500  jp = indxp( j )
501  dlamda( j ) = d( jp )
502  perm( j ) = indxq( indx( jp ) )
503  CALL dcopy( qsiz, q( 1, perm( j ) ), 1, q2( 1, j ), 1 )
504  130 CONTINUE
505  END IF
506 *
507 * The deflated eigenvalues and their corresponding vectors go back
508 * into the last N - K slots of D and Q respectively.
509 *
510  IF( k.LT.n ) THEN
511  IF( icompq.EQ.0 ) THEN
512  CALL dcopy( n-k, dlamda( k+1 ), 1, d( k+1 ), 1 )
513  ELSE
514  CALL dcopy( n-k, dlamda( k+1 ), 1, d( k+1 ), 1 )
515  CALL dlacpy( 'A', qsiz, n-k, q2( 1, k+1 ), ldq2,
516  $ q( 1, k+1 ), ldq )
517  END IF
518  END IF
519 *
520  RETURN
521 *
522 * End of DLAED8
523 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine drot(N, DX, INCX, DY, INCY, C, S)
DROT
Definition: drot.f:53
subroutine dlamrg(N1, N2, A, DTRD1, DTRD2, INDEX)
DLAMRG creates a permutation list to merge the entries of two independently sorted sets into a single...
Definition: dlamrg.f:101
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
double precision function dlapy2(X, Y)
DLAPY2 returns sqrt(x2+y2).
Definition: dlapy2.f:65

Here is the call graph for this function:

Here is the caller graph for this function: