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

◆ dlaed3()

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

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

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

Purpose:
 DLAED3 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 DLAED4 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
          DLAED4.  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 DOUBLE PRECISION array, dimension (N)
          D(I) contains the updated eigenvalues for
          1 <= I <= K.
[out]Q
          Q is DOUBLE PRECISION 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 DOUBLE PRECISION
          The value of the parameter in the rank one update equation.
          RHO >= 0 required.
[in]DLAMBDA
          DLAMBDA is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DLAED2).
          The rows of the eigenvectors found by DLAED4 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 dlaed3.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 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*
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 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
real(wp) function dnrm2(n, x, incx)
DNRM2
Definition dnrm2.f90:89
Here is the call graph for this function:
Here is the caller graph for this function: