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