LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlaed2.f
Go to the documentation of this file.
1*> \brief \b DLAED2 used by DSTEDC. Merges eigenvalues and deflates secular equation. Used when the original matrix is tridiagonal.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DLAED2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaed2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaed2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaed2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DLAED2( K, N, N1, D, Q, LDQ, INDXQ, RHO, Z, DLAMBDA, W,
20* Q2, INDX, INDXC, INDXP, COLTYP, INFO )
21*
22* .. Scalar Arguments ..
23* INTEGER INFO, K, LDQ, N, N1
24* DOUBLE PRECISION RHO
25* ..
26* .. Array Arguments ..
27* INTEGER COLTYP( * ), INDX( * ), INDXC( * ), INDXP( * ),
28* $ INDXQ( * )
29* DOUBLE PRECISION D( * ), DLAMBDA( * ), Q( LDQ, * ), Q2( * ),
30* $ W( * ), Z( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> DLAED2 merges the two sets of eigenvalues together into a single
40*> sorted set. Then it tries to deflate the size of the problem.
41*> There are two ways in which deflation can occur: when two or more
42*> eigenvalues are close together or if there is a tiny entry in the
43*> Z vector. For each such occurrence the order of the related secular
44*> equation problem is reduced by one.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[out] K
51*> \verbatim
52*> K is INTEGER
53*> The number of non-deflated eigenvalues, and the order of the
54*> related secular equation. 0 <= K <=N.
55*> \endverbatim
56*>
57*> \param[in] N
58*> \verbatim
59*> N is INTEGER
60*> The dimension of the symmetric tridiagonal matrix. N >= 0.
61*> \endverbatim
62*>
63*> \param[in] N1
64*> \verbatim
65*> N1 is INTEGER
66*> The location of the last eigenvalue in the leading sub-matrix.
67*> min(1,N) <= N1 <= N/2.
68*> \endverbatim
69*>
70*> \param[in,out] D
71*> \verbatim
72*> D is DOUBLE PRECISION array, dimension (N)
73*> On entry, D contains the eigenvalues of the two submatrices to
74*> be combined.
75*> On exit, D contains the trailing (N-K) updated eigenvalues
76*> (those which were deflated) sorted into increasing order.
77*> \endverbatim
78*>
79*> \param[in,out] Q
80*> \verbatim
81*> Q is DOUBLE PRECISION array, dimension (LDQ, N)
82*> On entry, Q contains the eigenvectors of two submatrices in
83*> the two square blocks with corners at (1,1), (N1,N1)
84*> and (N1+1, N1+1), (N,N).
85*> On exit, Q contains the trailing (N-K) updated eigenvectors
86*> (those which were deflated) in its last N-K columns.
87*> \endverbatim
88*>
89*> \param[in] LDQ
90*> \verbatim
91*> LDQ is INTEGER
92*> The leading dimension of the array Q. LDQ >= max(1,N).
93*> \endverbatim
94*>
95*> \param[in,out] INDXQ
96*> \verbatim
97*> INDXQ is INTEGER array, dimension (N)
98*> The permutation which separately sorts the two sub-problems
99*> in D into ascending order. Note that elements in the second
100*> half of this permutation must first have N1 added to their
101*> values. Destroyed on exit.
102*> \endverbatim
103*>
104*> \param[in,out] RHO
105*> \verbatim
106*> RHO is DOUBLE PRECISION
107*> On entry, the off-diagonal element associated with the rank-1
108*> cut which originally split the two submatrices which are now
109*> being recombined.
110*> On exit, RHO has been modified to the value required by
111*> DLAED3.
112*> \endverbatim
113*>
114*> \param[in] Z
115*> \verbatim
116*> Z is DOUBLE PRECISION array, dimension (N)
117*> On entry, Z contains the updating vector (the last
118*> row of the first sub-eigenvector matrix and the first row of
119*> the second sub-eigenvector matrix).
120*> On exit, the contents of Z have been destroyed by the updating
121*> process.
122*> \endverbatim
123*>
124*> \param[out] DLAMBDA
125*> \verbatim
126*> DLAMBDA is DOUBLE PRECISION array, dimension (N)
127*> A copy of the first K eigenvalues which will be used by
128*> DLAED3 to form the secular equation.
129*> \endverbatim
130*>
131*> \param[out] W
132*> \verbatim
133*> W is DOUBLE PRECISION array, dimension (N)
134*> The first k values of the final deflation-altered z-vector
135*> which will be passed to DLAED3.
136*> \endverbatim
137*>
138*> \param[out] Q2
139*> \verbatim
140*> Q2 is DOUBLE PRECISION array, dimension (N1**2+(N-N1)**2)
141*> A copy of the first K eigenvectors which will be used by
142*> DLAED3 in a matrix multiply (DGEMM) to solve for the new
143*> eigenvectors.
144*> \endverbatim
145*>
146*> \param[out] INDX
147*> \verbatim
148*> INDX is INTEGER array, dimension (N)
149*> The permutation used to sort the contents of DLAMBDA into
150*> ascending order.
151*> \endverbatim
152*>
153*> \param[out] INDXC
154*> \verbatim
155*> INDXC is INTEGER array, dimension (N)
156*> The permutation used to arrange the columns of the deflated
157*> Q matrix into three groups: the first group contains non-zero
158*> elements only at and above N1, the second contains
159*> non-zero elements only below N1, and the third is dense.
160*> \endverbatim
161*>
162*> \param[out] INDXP
163*> \verbatim
164*> INDXP is INTEGER array, dimension (N)
165*> The permutation used to place deflated values of D at the end
166*> of the array. INDXP(1:K) points to the nondeflated D-values
167*> and INDXP(K+1:N) points to the deflated eigenvalues.
168*> \endverbatim
169*>
170*> \param[out] COLTYP
171*> \verbatim
172*> COLTYP is INTEGER array, dimension (N)
173*> During execution, a label which will indicate which of the
174*> following types a column in the Q2 matrix is:
175*> 1 : non-zero in the upper half only;
176*> 2 : dense;
177*> 3 : non-zero in the lower half only;
178*> 4 : deflated.
179*> On exit, COLTYP(i) is the number of columns of type i,
180*> for i=1 to 4 only.
181*> \endverbatim
182*>
183*> \param[out] INFO
184*> \verbatim
185*> INFO is INTEGER
186*> = 0: successful exit.
187*> < 0: if INFO = -i, the i-th argument had an illegal value.
188*> \endverbatim
189*
190* Authors:
191* ========
192*
193*> \author Univ. of Tennessee
194*> \author Univ. of California Berkeley
195*> \author Univ. of Colorado Denver
196*> \author NAG Ltd.
197*
198*> \ingroup laed2
199*
200*> \par Contributors:
201* ==================
202*>
203*> Jeff Rutter, Computer Science Division, University of California
204*> at Berkeley, USA \n
205*> Modified by Francoise Tisseur, University of Tennessee
206*>
207* =====================================================================
208 SUBROUTINE dlaed2( K, N, N1, D, Q, LDQ, INDXQ, RHO, Z, DLAMBDA,
209 $ W,
210 $ Q2, INDX, INDXC, INDXP, COLTYP, INFO )
211*
212* -- LAPACK computational routine --
213* -- LAPACK is a software package provided by Univ. of Tennessee, --
214* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
215*
216* .. Scalar Arguments ..
217 INTEGER INFO, K, LDQ, N, N1
218 DOUBLE PRECISION RHO
219* ..
220* .. Array Arguments ..
221 INTEGER COLTYP( * ), INDX( * ), INDXC( * ), INDXP( * ),
222 $ INDXQ( * )
223 DOUBLE PRECISION D( * ), DLAMBDA( * ), Q( LDQ, * ), Q2( * ),
224 $ W( * ), Z( * )
225* ..
226*
227* =====================================================================
228*
229* .. Parameters ..
230 DOUBLE PRECISION MONE, ZERO, ONE, TWO, EIGHT
231 PARAMETER ( MONE = -1.0d0, zero = 0.0d0, one = 1.0d0,
232 $ two = 2.0d0, eight = 8.0d0 )
233* ..
234* .. Local Arrays ..
235 INTEGER CTOT( 4 ), PSM( 4 )
236* ..
237* .. Local Scalars ..
238 INTEGER CT, I, IMAX, IQ1, IQ2, J, JMAX, JS, K2, N1P1,
239 $ N2, NJ, PJ
240 DOUBLE PRECISION C, EPS, S, T, TAU, TOL
241* ..
242* .. External Functions ..
243 INTEGER IDAMAX
244 DOUBLE PRECISION DLAMCH, DLAPY2
245 EXTERNAL idamax, dlamch, dlapy2
246* ..
247* .. External Subroutines ..
248 EXTERNAL dcopy, dlacpy, dlamrg, drot, dscal,
249 $ xerbla
250* ..
251* .. Intrinsic Functions ..
252 INTRINSIC abs, max, min, sqrt
253* ..
254* .. Executable Statements ..
255*
256* Test the input parameters.
257*
258 info = 0
259*
260 IF( n.LT.0 ) THEN
261 info = -2
262 ELSE IF( ldq.LT.max( 1, n ) ) THEN
263 info = -6
264 ELSE IF( min( 1, ( n / 2 ) ).GT.n1 .OR. ( n / 2 ).LT.n1 ) THEN
265 info = -3
266 END IF
267 IF( info.NE.0 ) THEN
268 CALL xerbla( 'DLAED2', -info )
269 RETURN
270 END IF
271*
272* Quick return if possible
273*
274 IF( n.EQ.0 )
275 $ RETURN
276*
277 n2 = n - n1
278 n1p1 = n1 + 1
279*
280 IF( rho.LT.zero ) THEN
281 CALL dscal( n2, mone, z( n1p1 ), 1 )
282 END IF
283*
284* Normalize z so that norm(z) = 1. Since z is the concatenation of
285* two normalized vectors, norm2(z) = sqrt(2).
286*
287 t = one / sqrt( two )
288 CALL dscal( n, t, z, 1 )
289*
290* RHO = ABS( norm(z)**2 * RHO )
291*
292 rho = abs( two*rho )
293*
294* Sort the eigenvalues into increasing order
295*
296 DO 10 i = n1p1, n
297 indxq( i ) = indxq( i ) + n1
298 10 CONTINUE
299*
300* re-integrate the deflated parts from the last pass
301*
302 DO 20 i = 1, n
303 dlambda( i ) = d( indxq( i ) )
304 20 CONTINUE
305 CALL dlamrg( n1, n2, dlambda, 1, 1, indxc )
306 DO 30 i = 1, n
307 indx( i ) = indxq( indxc( i ) )
308 30 CONTINUE
309*
310* Calculate the allowable deflation tolerance
311*
312 imax = idamax( n, z, 1 )
313 jmax = idamax( n, d, 1 )
314 eps = dlamch( 'Epsilon' )
315 tol = eight*eps*max( abs( d( jmax ) ), abs( z( imax ) ) )
316*
317* If the rank-1 modifier is small enough, no more needs to be done
318* except to reorganize Q so that its columns correspond with the
319* elements in D.
320*
321 IF( rho*abs( z( imax ) ).LE.tol ) THEN
322 k = 0
323 iq2 = 1
324 DO 40 j = 1, n
325 i = indx( j )
326 CALL dcopy( n, q( 1, i ), 1, q2( iq2 ), 1 )
327 dlambda( j ) = d( i )
328 iq2 = iq2 + n
329 40 CONTINUE
330 CALL dlacpy( 'A', n, n, q2, n, q, ldq )
331 CALL dcopy( n, dlambda, 1, d, 1 )
332 GO TO 190
333 END IF
334*
335* If there are multiple eigenvalues then the problem deflates. Here
336* the number of equal eigenvalues are found. As each equal
337* eigenvalue is found, an elementary reflector is computed to rotate
338* the corresponding eigensubspace so that the corresponding
339* components of Z are zero in this new basis.
340*
341 DO 50 i = 1, n1
342 coltyp( i ) = 1
343 50 CONTINUE
344 DO 60 i = n1p1, n
345 coltyp( i ) = 3
346 60 CONTINUE
347*
348*
349 k = 0
350 k2 = n + 1
351 DO 70 j = 1, n
352 nj = indx( j )
353 IF( rho*abs( z( nj ) ).LE.tol ) THEN
354*
355* Deflate due to small z component.
356*
357 k2 = k2 - 1
358 coltyp( nj ) = 4
359 indxp( k2 ) = nj
360 IF( j.EQ.n )
361 $ GO TO 100
362 ELSE
363 pj = nj
364 GO TO 80
365 END IF
366 70 CONTINUE
367 80 CONTINUE
368 j = j + 1
369 nj = indx( j )
370 IF( j.GT.n )
371 $ GO TO 100
372 IF( rho*abs( z( nj ) ).LE.tol ) THEN
373*
374* Deflate due to small z component.
375*
376 k2 = k2 - 1
377 coltyp( nj ) = 4
378 indxp( k2 ) = nj
379 ELSE
380*
381* Check if eigenvalues are close enough to allow deflation.
382*
383 s = z( pj )
384 c = z( nj )
385*
386* Find sqrt(a**2+b**2) without overflow or
387* destructive underflow.
388*
389 tau = dlapy2( c, s )
390 t = d( nj ) - d( pj )
391 c = c / tau
392 s = -s / tau
393 IF( abs( t*c*s ).LE.tol ) THEN
394*
395* Deflation is possible.
396*
397 z( nj ) = tau
398 z( pj ) = zero
399 IF( coltyp( nj ).NE.coltyp( pj ) )
400 $ coltyp( nj ) = 2
401 coltyp( pj ) = 4
402 CALL drot( n, q( 1, pj ), 1, q( 1, nj ), 1, c, s )
403 t = d( pj )*c**2 + d( nj )*s**2
404 d( nj ) = d( pj )*s**2 + d( nj )*c**2
405 d( pj ) = t
406 k2 = k2 - 1
407 i = 1
408 90 CONTINUE
409 IF( k2+i.LE.n ) THEN
410 IF( d( pj ).LT.d( indxp( k2+i ) ) ) THEN
411 indxp( k2+i-1 ) = indxp( k2+i )
412 indxp( k2+i ) = pj
413 i = i + 1
414 GO TO 90
415 ELSE
416 indxp( k2+i-1 ) = pj
417 END IF
418 ELSE
419 indxp( k2+i-1 ) = pj
420 END IF
421 pj = nj
422 ELSE
423 k = k + 1
424 dlambda( k ) = d( pj )
425 w( k ) = z( pj )
426 indxp( k ) = pj
427 pj = nj
428 END IF
429 END IF
430 GO TO 80
431 100 CONTINUE
432*
433* Record the last eigenvalue.
434*
435 k = k + 1
436 dlambda( k ) = d( pj )
437 w( k ) = z( pj )
438 indxp( k ) = pj
439*
440* Count up the total number of the various types of columns, then
441* form a permutation which positions the four column types into
442* four uniform groups (although one or more of these groups may be
443* empty).
444*
445 DO 110 j = 1, 4
446 ctot( j ) = 0
447 110 CONTINUE
448 DO 120 j = 1, n
449 ct = coltyp( j )
450 ctot( ct ) = ctot( ct ) + 1
451 120 CONTINUE
452*
453* PSM(*) = Position in SubMatrix (of types 1 through 4)
454*
455 psm( 1 ) = 1
456 psm( 2 ) = 1 + ctot( 1 )
457 psm( 3 ) = psm( 2 ) + ctot( 2 )
458 psm( 4 ) = psm( 3 ) + ctot( 3 )
459 k = n - ctot( 4 )
460*
461* Fill out the INDXC array so that the permutation which it induces
462* will place all type-1 columns first, all type-2 columns next,
463* then all type-3's, and finally all type-4's.
464*
465 DO 130 j = 1, n
466 js = indxp( j )
467 ct = coltyp( js )
468 indx( psm( ct ) ) = js
469 indxc( psm( ct ) ) = j
470 psm( ct ) = psm( ct ) + 1
471 130 CONTINUE
472*
473* Sort the eigenvalues and corresponding eigenvectors into DLAMBDA
474* and Q2 respectively. The eigenvalues/vectors which were not
475* deflated go into the first K slots of DLAMBDA and Q2 respectively,
476* while those which were deflated go into the last N - K slots.
477*
478 i = 1
479 iq1 = 1
480 iq2 = 1 + ( ctot( 1 )+ctot( 2 ) )*n1
481 DO 140 j = 1, ctot( 1 )
482 js = indx( i )
483 CALL dcopy( n1, q( 1, js ), 1, q2( iq1 ), 1 )
484 z( i ) = d( js )
485 i = i + 1
486 iq1 = iq1 + n1
487 140 CONTINUE
488*
489 DO 150 j = 1, ctot( 2 )
490 js = indx( i )
491 CALL dcopy( n1, q( 1, js ), 1, q2( iq1 ), 1 )
492 CALL dcopy( n2, q( n1+1, js ), 1, q2( iq2 ), 1 )
493 z( i ) = d( js )
494 i = i + 1
495 iq1 = iq1 + n1
496 iq2 = iq2 + n2
497 150 CONTINUE
498*
499 DO 160 j = 1, ctot( 3 )
500 js = indx( i )
501 CALL dcopy( n2, q( n1+1, js ), 1, q2( iq2 ), 1 )
502 z( i ) = d( js )
503 i = i + 1
504 iq2 = iq2 + n2
505 160 CONTINUE
506*
507 iq1 = iq2
508 DO 170 j = 1, ctot( 4 )
509 js = indx( i )
510 CALL dcopy( n, q( 1, js ), 1, q2( iq2 ), 1 )
511 iq2 = iq2 + n
512 z( i ) = d( js )
513 i = i + 1
514 170 CONTINUE
515*
516* The deflated eigenvalues and their corresponding vectors go back
517* into the last N - K slots of D and Q respectively.
518*
519 IF( k.LT.n ) THEN
520 CALL dlacpy( 'A', n, ctot( 4 ), q2( iq1 ), n,
521 $ q( 1, k+1 ), ldq )
522 CALL dcopy( n-k, z( k+1 ), 1, d( k+1 ), 1 )
523 END IF
524*
525* Copy CTOT into COLTYP for referencing in DLAED3.
526*
527 DO 180 j = 1, 4
528 coltyp( j ) = ctot( j )
529 180 CONTINUE
530*
531 190 CONTINUE
532 RETURN
533*
534* End of DLAED2
535*
536 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
subroutine dlaed2(k, n, n1, d, q, ldq, indxq, rho, z, dlambda, w, q2, indx, indxc, indxp, coltyp, info)
DLAED2 used by DSTEDC. Merges eigenvalues and deflates secular equation. Used when the original matri...
Definition dlaed2.f:211
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:97
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