LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlaed9 ( integer  K,
integer  KSTART,
integer  KSTOP,
integer  N,
double precision, dimension( * )  D,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
double precision  RHO,
double precision, dimension( * )  DLAMDA,
double precision, dimension( * )  W,
double precision, dimension( lds, * )  S,
integer  LDS,
integer  INFO 
)

DLAED9 used by sstedc. Finds the roots of the secular equation and updates the eigenvectors. Used when the original matrix is dense.

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

Purpose:
 DLAED9 finds the roots of the secular equation, as defined by the
 values in D, Z, and RHO, between KSTART and KSTOP.  It makes the
 appropriate calls to DLAED4 and then stores the new matrix of
 eigenvectors for use in calculating the next level of Z vectors.
Parameters
[in]K
          K is INTEGER
          The number of terms in the rational function to be solved by
          DLAED4.  K >= 0.
[in]KSTART
          KSTART is INTEGER
[in]KSTOP
          KSTOP is INTEGER
          The updated eigenvalues Lambda(I), KSTART <= I <= KSTOP
          are to be computed.  1 <= KSTART <= KSTOP <= K.
[in]N
          N is INTEGER
          The number of rows and columns in the Q matrix.
          N >= K (delation may result in N > K).
[out]D
          D is DOUBLE PRECISION array, dimension (N)
          D(I) contains the updated eigenvalues
          for KSTART <= I <= KSTOP.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDQ,N)
[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]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.
[in]W
          W is DOUBLE PRECISION array, dimension (K)
          The first K elements of this array contain the components
          of the deflation-adjusted updating vector.
[out]S
          S is DOUBLE PRECISION array, dimension (LDS, K)
          Will contain the eigenvectors of the repaired matrix which
          will be stored for subsequent Z vector calculation and
          multiplied by the previously accumulated eigenvectors
          to update the system.
[in]LDS
          LDS is INTEGER
          The leading dimension of S.  LDS >= max( 1, K ).
[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.
Date
September 2012
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Definition at line 158 of file dlaed9.f.

158 *
159 * -- LAPACK computational routine (version 3.4.2) --
160 * -- LAPACK is a software package provided by Univ. of Tennessee, --
161 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
162 * September 2012
163 *
164 * .. Scalar Arguments ..
165  INTEGER info, k, kstart, kstop, ldq, lds, n
166  DOUBLE PRECISION rho
167 * ..
168 * .. Array Arguments ..
169  DOUBLE PRECISION d( * ), dlamda( * ), q( ldq, * ), s( lds, * ),
170  $ w( * )
171 * ..
172 *
173 * =====================================================================
174 *
175 * .. Local Scalars ..
176  INTEGER i, j
177  DOUBLE PRECISION temp
178 * ..
179 * .. External Functions ..
180  DOUBLE PRECISION dlamc3, dnrm2
181  EXTERNAL dlamc3, dnrm2
182 * ..
183 * .. External Subroutines ..
184  EXTERNAL dcopy, dlaed4, xerbla
185 * ..
186 * .. Intrinsic Functions ..
187  INTRINSIC max, sign, sqrt
188 * ..
189 * .. Executable Statements ..
190 *
191 * Test the input parameters.
192 *
193  info = 0
194 *
195  IF( k.LT.0 ) THEN
196  info = -1
197  ELSE IF( kstart.LT.1 .OR. kstart.GT.max( 1, k ) ) THEN
198  info = -2
199  ELSE IF( max( 1, kstop ).LT.kstart .OR. kstop.GT.max( 1, k ) )
200  $ THEN
201  info = -3
202  ELSE IF( n.LT.k ) THEN
203  info = -4
204  ELSE IF( ldq.LT.max( 1, k ) ) THEN
205  info = -7
206  ELSE IF( lds.LT.max( 1, k ) ) THEN
207  info = -12
208  END IF
209  IF( info.NE.0 ) THEN
210  CALL xerbla( 'DLAED9', -info )
211  RETURN
212  END IF
213 *
214 * Quick return if possible
215 *
216  IF( k.EQ.0 )
217  $ RETURN
218 *
219 * Modify values DLAMDA(i) to make sure all DLAMDA(i)-DLAMDA(j) can
220 * be computed with high relative accuracy (barring over/underflow).
221 * This is a problem on machines without a guard digit in
222 * add/subtract (Cray XMP, Cray YMP, Cray C 90 and Cray 2).
223 * The following code replaces DLAMDA(I) by 2*DLAMDA(I)-DLAMDA(I),
224 * which on any of these machines zeros out the bottommost
225 * bit of DLAMDA(I) if it is 1; this makes the subsequent
226 * subtractions DLAMDA(I)-DLAMDA(J) unproblematic when cancellation
227 * occurs. On binary machines with a guard digit (almost all
228 * machines) it does not change DLAMDA(I) at all. On hexadecimal
229 * and decimal machines with a guard digit, it slightly
230 * changes the bottommost bits of DLAMDA(I). It does not account
231 * for hexadecimal or decimal machines without guard digits
232 * (we know of none). We use a subroutine call to compute
233 * 2*DLAMBDA(I) to prevent optimizing compilers from eliminating
234 * this code.
235 *
236  DO 10 i = 1, n
237  dlamda( i ) = dlamc3( dlamda( i ), dlamda( i ) ) - dlamda( i )
238  10 CONTINUE
239 *
240  DO 20 j = kstart, kstop
241  CALL dlaed4( k, j, dlamda, w, q( 1, j ), rho, d( j ), info )
242 *
243 * If the zero finder fails, the computation is terminated.
244 *
245  IF( info.NE.0 )
246  $ GO TO 120
247  20 CONTINUE
248 *
249  IF( k.EQ.1 .OR. k.EQ.2 ) THEN
250  DO 40 i = 1, k
251  DO 30 j = 1, k
252  s( j, i ) = q( j, i )
253  30 CONTINUE
254  40 CONTINUE
255  GO TO 120
256  END IF
257 *
258 * Compute updated W.
259 *
260  CALL dcopy( k, w, 1, s, 1 )
261 *
262 * Initialize W(I) = Q(I,I)
263 *
264  CALL dcopy( k, q, ldq+1, w, 1 )
265  DO 70 j = 1, k
266  DO 50 i = 1, j - 1
267  w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
268  50 CONTINUE
269  DO 60 i = j + 1, k
270  w( i ) = w( i )*( q( i, j ) / ( dlamda( i )-dlamda( j ) ) )
271  60 CONTINUE
272  70 CONTINUE
273  DO 80 i = 1, k
274  w( i ) = sign( sqrt( -w( i ) ), s( i, 1 ) )
275  80 CONTINUE
276 *
277 * Compute eigenvectors of the modified rank-1 modification.
278 *
279  DO 110 j = 1, k
280  DO 90 i = 1, k
281  q( i, j ) = w( i ) / q( i, j )
282  90 CONTINUE
283  temp = dnrm2( k, q( 1, j ), 1 )
284  DO 100 i = 1, k
285  s( i, j ) = q( i, j ) / temp
286  100 CONTINUE
287  110 CONTINUE
288 *
289  120 CONTINUE
290  RETURN
291 *
292 * End of DLAED9
293 *
double precision function dlamc3(A, B)
DLAMC3
Definition: dlamch.f:169
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
double precision function dnrm2(N, X, INCX)
DNRM2
Definition: dnrm2.f:56
subroutine dlaed4(N, I, D, Z, DELTA, RHO, DLAM, INFO)
DLAED4 used by sstedc. Finds a single root of the secular equation.
Definition: dlaed4.f:147

Here is the call graph for this function:

Here is the caller graph for this function: