LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlaed3.f
Go to the documentation of this file.
1*> \brief \b DLAED3 used by DSTEDC. Finds the roots of the secular equation and updates the eigenvectors. 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*> \htmlonly
9*> Download DLAED3 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaed3.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaed3.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaed3.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DLAED3( K, N, N1, D, Q, LDQ, RHO, DLAMBDA, Q2, INDX,
22* CTOT, W, S, INFO )
23*
24* .. Scalar Arguments ..
25* INTEGER INFO, K, LDQ, N, N1
26* DOUBLE PRECISION RHO
27* ..
28* .. Array Arguments ..
29* INTEGER CTOT( * ), INDX( * )
30* DOUBLE PRECISION D( * ), DLAMBDA( * ), Q( LDQ, * ), Q2( * ),
31* $ S( * ), W( * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> DLAED3 finds the roots of the secular equation, as defined by the
41*> values in D, W, and RHO, between 1 and K. It makes the
42*> appropriate calls to DLAED4 and then updates the eigenvectors by
43*> multiplying the matrix of eigenvectors of the pair of eigensystems
44*> being combined by the matrix of eigenvectors of the K-by-K system
45*> which is solved here.
46*>
47*> \endverbatim
48*
49* Arguments:
50* ==========
51*
52*> \param[in] K
53*> \verbatim
54*> K is INTEGER
55*> The number of terms in the rational function to be solved by
56*> DLAED4. K >= 0.
57*> \endverbatim
58*>
59*> \param[in] N
60*> \verbatim
61*> N is INTEGER
62*> The number of rows and columns in the Q matrix.
63*> N >= K (deflation may result in N>K).
64*> \endverbatim
65*>
66*> \param[in] N1
67*> \verbatim
68*> N1 is INTEGER
69*> The location of the last eigenvalue in the leading submatrix.
70*> min(1,N) <= N1 <= N/2.
71*> \endverbatim
72*>
73*> \param[out] D
74*> \verbatim
75*> D is DOUBLE PRECISION array, dimension (N)
76*> D(I) contains the updated eigenvalues for
77*> 1 <= I <= K.
78*> \endverbatim
79*>
80*> \param[out] Q
81*> \verbatim
82*> Q is DOUBLE PRECISION array, dimension (LDQ,N)
83*> Initially the first K columns are used as workspace.
84*> On output the columns 1 to K contain
85*> the updated eigenvectors.
86*> \endverbatim
87*>
88*> \param[in] LDQ
89*> \verbatim
90*> LDQ is INTEGER
91*> The leading dimension of the array Q. LDQ >= max(1,N).
92*> \endverbatim
93*>
94*> \param[in] RHO
95*> \verbatim
96*> RHO is DOUBLE PRECISION
97*> The value of the parameter in the rank one update equation.
98*> RHO >= 0 required.
99*> \endverbatim
100*>
101*> \param[in] DLAMBDA
102*> \verbatim
103*> DLAMBDA is DOUBLE PRECISION array, dimension (K)
104*> The first K elements of this array contain the old roots
105*> of the deflated updating problem. These are the poles
106*> of the secular equation.
107*> \endverbatim
108*>
109*> \param[in] Q2
110*> \verbatim
111*> Q2 is DOUBLE PRECISION array, dimension (LDQ2*N)
112*> The first K columns of this matrix contain the non-deflated
113*> eigenvectors for the split problem.
114*> \endverbatim
115*>
116*> \param[in] INDX
117*> \verbatim
118*> INDX is INTEGER array, dimension (N)
119*> The permutation used to arrange the columns of the deflated
120*> Q matrix into three groups (see DLAED2).
121*> The rows of the eigenvectors found by DLAED4 must be likewise
122*> permuted before the matrix multiply can take place.
123*> \endverbatim
124*>
125*> \param[in] CTOT
126*> \verbatim
127*> CTOT is INTEGER array, dimension (4)
128*> A count of the total number of the various types of columns
129*> in Q, as described in INDX. The fourth column type is any
130*> column which has been deflated.
131*> \endverbatim
132*>
133*> \param[in,out] W
134*> \verbatim
135*> W is DOUBLE PRECISION array, dimension (K)
136*> The first K elements of this array contain the components
137*> of the deflation-adjusted updating vector. Destroyed on
138*> output.
139*> \endverbatim
140*>
141*> \param[out] S
142*> \verbatim
143*> S is DOUBLE PRECISION array, dimension (N1 + 1)*K
144*> Will contain the eigenvectors of the repaired matrix which
145*> will be multiplied by the previously accumulated eigenvectors
146*> to update the system.
147*> \endverbatim
148*>
149*> \param[out] INFO
150*> \verbatim
151*> INFO is INTEGER
152*> = 0: successful exit.
153*> < 0: if INFO = -i, the i-th argument had an illegal value.
154*> > 0: if INFO = 1, an eigenvalue did not converge
155*> \endverbatim
156*
157* Authors:
158* ========
159*
160*> \author Univ. of Tennessee
161*> \author Univ. of California Berkeley
162*> \author Univ. of Colorado Denver
163*> \author NAG Ltd.
164*
165*> \ingroup laed3
166*
167*> \par Contributors:
168* ==================
169*>
170*> Jeff Rutter, Computer Science Division, University of California
171*> at Berkeley, USA \n
172*> Modified by Francoise Tisseur, University of Tennessee
173*>
174* =====================================================================
175 SUBROUTINE dlaed3( K, N, N1, D, Q, LDQ, RHO, DLAMBDA, Q2, INDX,
176 $ CTOT, W, S, INFO )
177*
178* -- LAPACK computational routine --
179* -- LAPACK is a software package provided by Univ. of Tennessee, --
180* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182* .. Scalar Arguments ..
183 INTEGER INFO, K, LDQ, N, N1
184 DOUBLE PRECISION RHO
185* ..
186* .. Array Arguments ..
187 INTEGER CTOT( * ), INDX( * )
188 DOUBLE PRECISION D( * ), DLAMBDA( * ), Q( LDQ, * ), Q2( * ),
189 $ s( * ), w( * )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 DOUBLE PRECISION ONE, ZERO
196 parameter( one = 1.0d0, zero = 0.0d0 )
197* ..
198* .. Local Scalars ..
199 INTEGER I, II, IQ2, J, N12, N2, N23
200 DOUBLE PRECISION TEMP
201* ..
202* .. External Functions ..
203 DOUBLE PRECISION DNRM2
204 EXTERNAL dnrm2
205* ..
206* .. External Subroutines ..
207 EXTERNAL dcopy, dgemm, dlacpy, dlaed4, dlaset, xerbla
208* ..
209* .. Intrinsic Functions ..
210 INTRINSIC max, sign, sqrt
211* ..
212* .. Executable Statements ..
213*
214* Test the input parameters.
215*
216 info = 0
217*
218 IF( k.LT.0 ) THEN
219 info = -1
220 ELSE IF( n.LT.k ) THEN
221 info = -2
222 ELSE IF( ldq.LT.max( 1, n ) ) THEN
223 info = -6
224 END IF
225 IF( info.NE.0 ) THEN
226 CALL xerbla( 'DLAED3', -info )
227 RETURN
228 END IF
229*
230* Quick return if possible
231*
232 IF( k.EQ.0 )
233 $ RETURN
234*
235*
236 DO 20 j = 1, k
237 CALL dlaed4( k, j, dlambda, w, q( 1, j ), rho, d( j ), info )
238*
239* If the zero finder fails, the computation is terminated.
240*
241 IF( info.NE.0 )
242 $ GO TO 120
243 20 CONTINUE
244*
245 IF( k.EQ.1 )
246 $ GO TO 110
247 IF( k.EQ.2 ) THEN
248 DO 30 j = 1, k
249 w( 1 ) = q( 1, j )
250 w( 2 ) = q( 2, j )
251 ii = indx( 1 )
252 q( 1, j ) = w( ii )
253 ii = indx( 2 )
254 q( 2, j ) = w( ii )
255 30 CONTINUE
256 GO TO 110
257 END IF
258*
259* Compute updated W.
260*
261 CALL dcopy( k, w, 1, s, 1 )
262*
263* Initialize W(I) = Q(I,I)
264*
265 CALL dcopy( k, q, ldq+1, w, 1 )
266 DO 60 j = 1, k
267 DO 40 i = 1, j - 1
268 w( i ) = w( i )*( q( i, j )/( dlambda( i )-dlambda( j ) ) )
269 40 CONTINUE
270 DO 50 i = j + 1, k
271 w( i ) = w( i )*( q( i, j )/( dlambda( i )-dlambda( j ) ) )
272 50 CONTINUE
273 60 CONTINUE
274 DO 70 i = 1, k
275 w( i ) = sign( sqrt( -w( i ) ), s( i ) )
276 70 CONTINUE
277*
278* Compute eigenvectors of the modified rank-1 modification.
279*
280 DO 100 j = 1, k
281 DO 80 i = 1, k
282 s( i ) = w( i ) / q( i, j )
283 80 CONTINUE
284 temp = dnrm2( k, s, 1 )
285 DO 90 i = 1, k
286 ii = indx( i )
287 q( i, j ) = s( ii ) / temp
288 90 CONTINUE
289 100 CONTINUE
290*
291* Compute the updated eigenvectors.
292*
293 110 CONTINUE
294*
295 n2 = n - n1
296 n12 = ctot( 1 ) + ctot( 2 )
297 n23 = ctot( 2 ) + ctot( 3 )
298*
299 CALL dlacpy( 'A', n23, k, q( ctot( 1 )+1, 1 ), ldq, s, n23 )
300 iq2 = n1*n12 + 1
301 IF( n23.NE.0 ) THEN
302 CALL dgemm( 'N', 'N', n2, k, n23, one, q2( iq2 ), n2, s, n23,
303 $ zero, q( n1+1, 1 ), ldq )
304 ELSE
305 CALL dlaset( 'A', n2, k, zero, zero, q( n1+1, 1 ), ldq )
306 END IF
307*
308 CALL dlacpy( 'A', n12, k, q, ldq, s, n12 )
309 IF( n12.NE.0 ) THEN
310 CALL dgemm( 'N', 'N', n1, k, n12, one, q2, n1, s, n12, zero, q,
311 $ ldq )
312 ELSE
313 CALL dlaset( 'A', n1, k, zero, zero, q( 1, 1 ), ldq )
314 END IF
315*
316*
317 120 CONTINUE
318 RETURN
319*
320* End of DLAED3
321*
322 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
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
subroutine dlaed3(k, n, n1, d, q, ldq, rho, dlambda, q2, indx, ctot, w, s, info)
DLAED3 used by DSTEDC. Finds the roots of the secular equation and updates the eigenvectors....
Definition dlaed3.f:177
subroutine dlaed4(n, i, d, z, delta, rho, dlam, info)
DLAED4 used by DSTEDC. Finds a single root of the secular equation.
Definition dlaed4.f:145
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:110