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

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

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

Purpose:
 CLAED8 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
[out]K
          K is INTEGER
         Contains the number of non-deflated eigenvalues.
         This is 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 unitary matrix used to reduce
         the dense or band matrix to tridiagonal form.
         QSIZ >= N if ICOMPQ = 1.
[in,out]Q
          Q is COMPLEX array, dimension (LDQ,N)
         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,out]D
          D is REAL array, dimension (N)
         On entry, D contains the eigenvalues of the two submatrices to
         be combined.  On exit, D contains the trailing (N-K) updated
         eigenvalues (those which were deflated) sorted into increasing
         order.
[in,out]RHO
          RHO is REAL
         Contains the off diagonal element associated with the rank-1
         cut which originally split the two submatrices which are now
         being recombined. RHO is modified during the computation to
         the value required by SLAED3.
[in]CUTPNT
          CUTPNT is INTEGER
         Contains the location of the last eigenvalue in the leading
         sub-matrix.  MIN(1,N) <= CUTPNT <= N.
[in]Z
          Z is REAL array, dimension (N)
         On input this vector contains the updating vector (the last
         row of the first sub-eigenvector matrix and the first row of
         the second sub-eigenvector matrix).  The contents of Z are
         destroyed during the updating process.
[out]DLAMDA
          DLAMDA is REAL array, dimension (N)
         Contains a copy of the first K eigenvalues which will be used
         by SLAED3 to form the secular equation.
[out]Q2
          Q2 is COMPLEX array, dimension (LDQ2,N)
         If ICOMPQ = 0, Q2 is not referenced.  Otherwise,
         Contains a copy of the first K eigenvectors which will be used
         by SLAED7 in a matrix multiply (SGEMM) 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 REAL array, dimension (N)
         This will hold the first k values of the final
         deflation-altered z-vector and will be passed to SLAED3.
[out]INDXP
          INDXP is INTEGER array, dimension (N)
         This will contain the permutation used to place deflated
         values of D at the end of the array. On output 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)
         This will contain the permutation used to sort the contents of
         D into ascending order.
[in]INDXQ
          INDXQ is INTEGER array, dimension (N)
         This contains 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.
[out]PERM
          PERM is INTEGER array, dimension (N)
         Contains the permutations (from deflation and sorting) to be
         applied to each eigenblock.
[out]GIVPTR
          GIVPTR is INTEGER
         Contains 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 REAL array, dimension (2, N)
         Each number indicates the S value to be used in the
         corresponding Givens rotation.
[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

Definition at line 230 of file claed8.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: