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

◆ slaed3()

subroutine slaed3 ( integer  k,
integer  n,
integer  n1,
real, dimension( * )  d,
real, dimension( ldq, * )  q,
integer  ldq,
real  rho,
real, dimension( * )  dlambda,
real, dimension( * )  q2,
integer, dimension( * )  indx,
integer, dimension( * )  ctot,
real, dimension( * )  w,
real, dimension( * )  s,
integer  info 
)

SLAED3 used by SSTEDC. Finds the roots of the secular equation and updates the eigenvectors. Used when the original matrix is tridiagonal.

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

Purpose:
 SLAED3 finds the roots of the secular equation, as defined by the
 values in D, W, and RHO, between 1 and K.  It makes the
 appropriate calls to SLAED4 and then updates the eigenvectors by
 multiplying the matrix of eigenvectors of the pair of eigensystems
 being combined by the matrix of eigenvectors of the K-by-K system
 which is solved here.
Parameters
[in]K
          K is INTEGER
          The number of terms in the rational function to be solved by
          SLAED4.  K >= 0.
[in]N
          N is INTEGER
          The number of rows and columns in the Q matrix.
          N >= K (deflation may result in N>K).
[in]N1
          N1 is INTEGER
          The location of the last eigenvalue in the leading submatrix.
          min(1,N) <= N1 <= N/2.
[out]D
          D is REAL array, dimension (N)
          D(I) contains the updated eigenvalues for
          1 <= I <= K.
[out]Q
          Q is REAL array, dimension (LDQ,N)
          Initially the first K columns are used as workspace.
          On output the columns 1 to K contain
          the updated eigenvectors.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  LDQ >= max(1,N).
[in]RHO
          RHO is REAL
          The value of the parameter in the rank one update equation.
          RHO >= 0 required.
[in]DLAMBDA
          DLAMBDA is REAL array, dimension (K)
          The first K elements of this array contain the old roots
          of the deflated updating problem.  These are the poles
          of the secular equation.
[in]Q2
          Q2 is REAL array, dimension (LDQ2*N)
          The first K columns of this matrix contain the non-deflated
          eigenvectors for the split problem.
[in]INDX
          INDX is INTEGER array, dimension (N)
          The permutation used to arrange the columns of the deflated
          Q matrix into three groups (see SLAED2).
          The rows of the eigenvectors found by SLAED4 must be likewise
          permuted before the matrix multiply can take place.
[in]CTOT
          CTOT is INTEGER array, dimension (4)
          A count of the total number of the various types of columns
          in Q, as described in INDX.  The fourth column type is any
          column which has been deflated.
[in,out]W
          W is REAL array, dimension (K)
          The first K elements of this array contain the components
          of the deflation-adjusted updating vector. Destroyed on
          output.
[out]S
          S is REAL array, dimension (N1 + 1)*K
          Will contain the eigenvectors of the repaired matrix which
          will be multiplied by the previously accumulated eigenvectors
          to update the system.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = 1, an eigenvalue did not converge
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
Modified by Francoise Tisseur, University of Tennessee

Definition at line 175 of file slaed3.f.

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 REAL RHO
185* ..
186* .. Array Arguments ..
187 INTEGER CTOT( * ), INDX( * )
188 REAL D( * ), DLAMBDA( * ), Q( LDQ, * ), Q2( * ),
189 $ S( * ), W( * )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 REAL ONE, ZERO
196 parameter( one = 1.0e0, zero = 0.0e0 )
197* ..
198* .. Local Scalars ..
199 INTEGER I, II, IQ2, J, N12, N2, N23
200 REAL TEMP
201* ..
202* .. External Functions ..
203 REAL SNRM2
204 EXTERNAL snrm2
205* ..
206* .. External Subroutines ..
207 EXTERNAL scopy, sgemm, slacpy, slaed4, slaset, 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( 'SLAED3', -info )
227 RETURN
228 END IF
229*
230* Quick return if possible
231*
232 IF( k.EQ.0 )
233 $ RETURN
234*
235 DO 20 j = 1, k
236 CALL slaed4( k, j, dlambda, w, q( 1, j ), rho, d( j ), 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, n23,
302 $ zero, q( n1+1, 1 ), ldq )
303 ELSE
304 CALL slaset( 'A', n2, k, zero, zero, q( n1+1, 1 ), ldq )
305 END IF
306*
307 CALL slacpy( 'A', n12, k, q, ldq, s, n12 )
308 IF( n12.NE.0 ) THEN
309 CALL sgemm( 'N', 'N', n1, k, n12, one, q2, n1, s, n12, zero, q,
310 $ ldq )
311 ELSE
312 CALL slaset( 'A', n1, k, zero, zero, q( 1, 1 ), ldq )
313 END IF
314*
315*
316 120 CONTINUE
317 RETURN
318*
319* End of SLAED3
320*
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:103
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:145
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:110
real(wp) function snrm2(n, x, incx)
SNRM2
Definition snrm2.f90:89
Here is the call graph for this function:
Here is the caller graph for this function: