LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dlaed8()

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( * )  dlambda,
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 DSTEDC. 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]DLAMBDA
          DLAMBDA 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.
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Definition at line 240 of file dlaed8.f.

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