LAPACK 3.11.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( * )  DLAMDA,
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.

 This code makes very mild assumptions about floating point
 arithmetic. It will work on machines with a guard digit in
 add/subtract, or on those binary machines without guard digits
 which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2.
 It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
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,out]DLAMDA
          DLAMDA 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. May be changed on output by
          having lowest order bit set to zero on Cray X-MP, Cray Y-MP,
          Cray-2, or Cray C-90, as described above.
[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 183 of file dlaed3.f.

185*
186* -- LAPACK computational routine --
187* -- LAPACK is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 INTEGER INFO, K, LDQ, N, N1
192 DOUBLE PRECISION RHO
193* ..
194* .. Array Arguments ..
195 INTEGER CTOT( * ), INDX( * )
196 DOUBLE PRECISION D( * ), DLAMDA( * ), Q( LDQ, * ), Q2( * ),
197 $ S( * ), W( * )
198* ..
199*
200* =====================================================================
201*
202* .. Parameters ..
203 DOUBLE PRECISION ONE, ZERO
204 parameter( one = 1.0d0, zero = 0.0d0 )
205* ..
206* .. Local Scalars ..
207 INTEGER I, II, IQ2, J, N12, N2, N23
208 DOUBLE PRECISION TEMP
209* ..
210* .. External Functions ..
211 DOUBLE PRECISION DLAMC3, DNRM2
212 EXTERNAL dlamc3, dnrm2
213* ..
214* .. External Subroutines ..
215 EXTERNAL dcopy, dgemm, dlacpy, dlaed4, dlaset, xerbla
216* ..
217* .. Intrinsic Functions ..
218 INTRINSIC max, sign, sqrt
219* ..
220* .. Executable Statements ..
221*
222* Test the input parameters.
223*
224 info = 0
225*
226 IF( k.LT.0 ) THEN
227 info = -1
228 ELSE IF( n.LT.k ) THEN
229 info = -2
230 ELSE IF( ldq.LT.max( 1, n ) ) THEN
231 info = -6
232 END IF
233 IF( info.NE.0 ) THEN
234 CALL xerbla( 'DLAED3', -info )
235 RETURN
236 END IF
237*
238* Quick return if possible
239*
240 IF( k.EQ.0 )
241 $ RETURN
242*
243* Modify values DLAMDA(i) to make sure all DLAMDA(i)-DLAMDA(j) can
244* be computed with high relative accuracy (barring over/underflow).
245* This is a problem on machines without a guard digit in
246* add/subtract (Cray XMP, Cray YMP, Cray C 90 and Cray 2).
247* The following code replaces DLAMDA(I) by 2*DLAMDA(I)-DLAMDA(I),
248* which on any of these machines zeros out the bottommost
249* bit of DLAMDA(I) if it is 1; this makes the subsequent
250* subtractions DLAMDA(I)-DLAMDA(J) unproblematic when cancellation
251* occurs. On binary machines with a guard digit (almost all
252* machines) it does not change DLAMDA(I) at all. On hexadecimal
253* and decimal machines with a guard digit, it slightly
254* changes the bottommost bits of DLAMDA(I). It does not account
255* for hexadecimal or decimal machines without guard digits
256* (we know of none). We use a subroutine call to compute
257* 2*DLAMBDA(I) to prevent optimizing compilers from eliminating
258* this code.
259*
260 DO 10 i = 1, k
261 dlamda( i ) = dlamc3( dlamda( i ), dlamda( i ) ) - dlamda( i )
262 10 CONTINUE
263*
264 DO 20 j = 1, k
265 CALL dlaed4( k, j, dlamda, w, q( 1, j ), rho, d( j ), info )
266*
267* If the zero finder fails, the computation is terminated.
268*
269 IF( info.NE.0 )
270 $ GO TO 120
271 20 CONTINUE
272*
273 IF( k.EQ.1 )
274 $ GO TO 110
275 IF( k.EQ.2 ) THEN
276 DO 30 j = 1, k
277 w( 1 ) = q( 1, j )
278 w( 2 ) = q( 2, j )
279 ii = indx( 1 )
280 q( 1, j ) = w( ii )
281 ii = indx( 2 )
282 q( 2, j ) = w( ii )
283 30 CONTINUE
284 GO TO 110
285 END IF
286*
287* Compute updated W.
288*
289 CALL dcopy( k, w, 1, s, 1 )
290*
291* Initialize W(I) = Q(I,I)
292*
293 CALL dcopy( k, q, ldq+1, w, 1 )
294 DO 60 j = 1, k
295 DO 40 i = 1, j - 1
296 w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
297 40 CONTINUE
298 DO 50 i = j + 1, k
299 w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
300 50 CONTINUE
301 60 CONTINUE
302 DO 70 i = 1, k
303 w( i ) = sign( sqrt( -w( i ) ), s( i ) )
304 70 CONTINUE
305*
306* Compute eigenvectors of the modified rank-1 modification.
307*
308 DO 100 j = 1, k
309 DO 80 i = 1, k
310 s( i ) = w( i ) / q( i, j )
311 80 CONTINUE
312 temp = dnrm2( k, s, 1 )
313 DO 90 i = 1, k
314 ii = indx( i )
315 q( i, j ) = s( ii ) / temp
316 90 CONTINUE
317 100 CONTINUE
318*
319* Compute the updated eigenvectors.
320*
321 110 CONTINUE
322*
323 n2 = n - n1
324 n12 = ctot( 1 ) + ctot( 2 )
325 n23 = ctot( 2 ) + ctot( 3 )
326*
327 CALL dlacpy( 'A', n23, k, q( ctot( 1 )+1, 1 ), ldq, s, n23 )
328 iq2 = n1*n12 + 1
329 IF( n23.NE.0 ) THEN
330 CALL dgemm( 'N', 'N', n2, k, n23, one, q2( iq2 ), n2, s, n23,
331 $ zero, q( n1+1, 1 ), ldq )
332 ELSE
333 CALL dlaset( 'A', n2, k, zero, zero, q( n1+1, 1 ), ldq )
334 END IF
335*
336 CALL dlacpy( 'A', n12, k, q, ldq, s, n12 )
337 IF( n12.NE.0 ) THEN
338 CALL dgemm( 'N', 'N', n1, k, n12, one, q2, n1, s, n12, zero, q,
339 $ ldq )
340 ELSE
341 CALL dlaset( 'A', n1, k, zero, zero, q( 1, 1 ), ldq )
342 END IF
343*
344*
345 120 CONTINUE
346 RETURN
347*
348* End of DLAED3
349*
double precision function dlamc3(A, B)
DLAMC3
Definition: dlamch.f:169
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 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
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
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 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:187
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: