LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sstein ( integer  N,
real, dimension( * )  D,
real, dimension( * )  E,
integer  M,
real, dimension( * )  W,
integer, dimension( * )  IBLOCK,
integer, dimension( * )  ISPLIT,
real, dimension( ldz, * )  Z,
integer  LDZ,
real, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

SSTEIN

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

Purpose:
 SSTEIN computes the eigenvectors of a real symmetric tridiagonal
 matrix T corresponding to specified eigenvalues, using inverse
 iteration.

 The maximum number of iterations allowed for each eigenvector is
 specified by an internal parameter MAXITS (currently set to 5).
Parameters
[in]N
          N is INTEGER
          The order of the matrix.  N >= 0.
[in]D
          D is REAL array, dimension (N)
          The n diagonal elements of the tridiagonal matrix T.
[in]E
          E is REAL array, dimension (N-1)
          The (n-1) subdiagonal elements of the tridiagonal matrix
          T, in elements 1 to N-1.
[in]M
          M is INTEGER
          The number of eigenvectors to be found.  0 <= M <= N.
[in]W
          W is REAL array, dimension (N)
          The first M elements of W contain the eigenvalues for
          which eigenvectors are to be computed.  The eigenvalues
          should be grouped by split-off block and ordered from
          smallest to largest within the block.  ( The output array
          W from SSTEBZ with ORDER = 'B' is expected here. )
[in]IBLOCK
          IBLOCK is INTEGER array, dimension (N)
          The submatrix indices associated with the corresponding
          eigenvalues in W; IBLOCK(i)=1 if eigenvalue W(i) belongs to
          the first submatrix from the top, =2 if W(i) belongs to
          the second submatrix, etc.  ( The output array IBLOCK
          from SSTEBZ is expected here. )
[in]ISPLIT
          ISPLIT is INTEGER array, dimension (N)
          The splitting points, at which T breaks up into submatrices.
          The first submatrix consists of rows/columns 1 to
          ISPLIT( 1 ), the second of rows/columns ISPLIT( 1 )+1
          through ISPLIT( 2 ), etc.
          ( The output array ISPLIT from SSTEBZ is expected here. )
[out]Z
          Z is REAL array, dimension (LDZ, M)
          The computed eigenvectors.  The eigenvector associated
          with the eigenvalue W(i) is stored in the i-th column of
          Z.  Any vector which fails to converge is set to its current
          iterate after MAXITS iterations.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= max(1,N).
[out]WORK
          WORK is REAL array, dimension (5*N)
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]IFAIL
          IFAIL is INTEGER array, dimension (M)
          On normal exit, all elements of IFAIL are zero.
          If one or more eigenvectors fail to converge after
          MAXITS iterations, then their indices are stored in
          array IFAIL.
[out]INFO
          INFO is INTEGER
          = 0: successful exit.
          < 0: if INFO = -i, the i-th argument had an illegal value
          > 0: if INFO = i, then i eigenvectors failed to converge
               in MAXITS iterations.  Their indices are stored in
               array IFAIL.
Internal Parameters:
  MAXITS  INTEGER, default = 5
          The maximum number of iterations performed.

  EXTRA   INTEGER, default = 2
          The number of iterations performed after norm growth
          criterion is satisfied, should be at least 1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 176 of file sstein.f.

176 *
177 * -- LAPACK computational routine (version 3.6.0) --
178 * -- LAPACK is a software package provided by Univ. of Tennessee, --
179 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
180 * November 2015
181 *
182 * .. Scalar Arguments ..
183  INTEGER info, ldz, m, n
184 * ..
185 * .. Array Arguments ..
186  INTEGER iblock( * ), ifail( * ), isplit( * ),
187  $ iwork( * )
188  REAL d( * ), e( * ), w( * ), work( * ), z( ldz, * )
189 * ..
190 *
191 * =====================================================================
192 *
193 * .. Parameters ..
194  REAL zero, one, ten, odm3, odm1
195  parameter ( zero = 0.0e+0, one = 1.0e+0, ten = 1.0e+1,
196  $ odm3 = 1.0e-3, odm1 = 1.0e-1 )
197  INTEGER maxits, extra
198  parameter ( maxits = 5, extra = 2 )
199 * ..
200 * .. Local Scalars ..
201  INTEGER b1, blksiz, bn, gpind, i, iinfo, indrv1,
202  $ indrv2, indrv3, indrv4, indrv5, its, j, j1,
203  $ jblk, jmax, nblk, nrmchk
204  REAL ctr, eps, eps1, nrm, onenrm, ortol, pertol,
205  $ scl, sep, stpcrt, tol, xj, xjm
206 * ..
207 * .. Local Arrays ..
208  INTEGER iseed( 4 )
209 * ..
210 * .. External Functions ..
211  INTEGER isamax
212  REAL sasum, sdot, slamch, snrm2
213  EXTERNAL isamax, sasum, sdot, slamch, snrm2
214 * ..
215 * .. External Subroutines ..
216  EXTERNAL saxpy, scopy, slagtf, slagts, slarnv, sscal,
217  $ xerbla
218 * ..
219 * .. Intrinsic Functions ..
220  INTRINSIC abs, max, sqrt
221 * ..
222 * .. Executable Statements ..
223 *
224 * Test the input parameters.
225 *
226  info = 0
227  DO 10 i = 1, m
228  ifail( i ) = 0
229  10 CONTINUE
230 *
231  IF( n.LT.0 ) THEN
232  info = -1
233  ELSE IF( m.LT.0 .OR. m.GT.n ) THEN
234  info = -4
235  ELSE IF( ldz.LT.max( 1, n ) ) THEN
236  info = -9
237  ELSE
238  DO 20 j = 2, m
239  IF( iblock( j ).LT.iblock( j-1 ) ) THEN
240  info = -6
241  GO TO 30
242  END IF
243  IF( iblock( j ).EQ.iblock( j-1 ) .AND. w( j ).LT.w( j-1 ) )
244  $ THEN
245  info = -5
246  GO TO 30
247  END IF
248  20 CONTINUE
249  30 CONTINUE
250  END IF
251 *
252  IF( info.NE.0 ) THEN
253  CALL xerbla( 'SSTEIN', -info )
254  RETURN
255  END IF
256 *
257 * Quick return if possible
258 *
259  IF( n.EQ.0 .OR. m.EQ.0 ) THEN
260  RETURN
261  ELSE IF( n.EQ.1 ) THEN
262  z( 1, 1 ) = one
263  RETURN
264  END IF
265 *
266 * Get machine constants.
267 *
268  eps = slamch( 'Precision' )
269 *
270 * Initialize seed for random number generator SLARNV.
271 *
272  DO 40 i = 1, 4
273  iseed( i ) = 1
274  40 CONTINUE
275 *
276 * Initialize pointers.
277 *
278  indrv1 = 0
279  indrv2 = indrv1 + n
280  indrv3 = indrv2 + n
281  indrv4 = indrv3 + n
282  indrv5 = indrv4 + n
283 *
284 * Compute eigenvectors of matrix blocks.
285 *
286  j1 = 1
287  DO 160 nblk = 1, iblock( m )
288 *
289 * Find starting and ending indices of block nblk.
290 *
291  IF( nblk.EQ.1 ) THEN
292  b1 = 1
293  ELSE
294  b1 = isplit( nblk-1 ) + 1
295  END IF
296  bn = isplit( nblk )
297  blksiz = bn - b1 + 1
298  IF( blksiz.EQ.1 )
299  $ GO TO 60
300  gpind = j1
301 *
302 * Compute reorthogonalization criterion and stopping criterion.
303 *
304  onenrm = abs( d( b1 ) ) + abs( e( b1 ) )
305  onenrm = max( onenrm, abs( d( bn ) )+abs( e( bn-1 ) ) )
306  DO 50 i = b1 + 1, bn - 1
307  onenrm = max( onenrm, abs( d( i ) )+abs( e( i-1 ) )+
308  $ abs( e( i ) ) )
309  50 CONTINUE
310  ortol = odm3*onenrm
311 *
312  stpcrt = sqrt( odm1 / blksiz )
313 *
314 * Loop through eigenvalues of block nblk.
315 *
316  60 CONTINUE
317  jblk = 0
318  DO 150 j = j1, m
319  IF( iblock( j ).NE.nblk ) THEN
320  j1 = j
321  GO TO 160
322  END IF
323  jblk = jblk + 1
324  xj = w( j )
325 *
326 * Skip all the work if the block size is one.
327 *
328  IF( blksiz.EQ.1 ) THEN
329  work( indrv1+1 ) = one
330  GO TO 120
331  END IF
332 *
333 * If eigenvalues j and j-1 are too close, add a relatively
334 * small perturbation.
335 *
336  IF( jblk.GT.1 ) THEN
337  eps1 = abs( eps*xj )
338  pertol = ten*eps1
339  sep = xj - xjm
340  IF( sep.LT.pertol )
341  $ xj = xjm + pertol
342  END IF
343 *
344  its = 0
345  nrmchk = 0
346 *
347 * Get random starting vector.
348 *
349  CALL slarnv( 2, iseed, blksiz, work( indrv1+1 ) )
350 *
351 * Copy the matrix T so it won't be destroyed in factorization.
352 *
353  CALL scopy( blksiz, d( b1 ), 1, work( indrv4+1 ), 1 )
354  CALL scopy( blksiz-1, e( b1 ), 1, work( indrv2+2 ), 1 )
355  CALL scopy( blksiz-1, e( b1 ), 1, work( indrv3+1 ), 1 )
356 *
357 * Compute LU factors with partial pivoting ( PT = LU )
358 *
359  tol = zero
360  CALL slagtf( blksiz, work( indrv4+1 ), xj, work( indrv2+2 ),
361  $ work( indrv3+1 ), tol, work( indrv5+1 ), iwork,
362  $ iinfo )
363 *
364 * Update iteration count.
365 *
366  70 CONTINUE
367  its = its + 1
368  IF( its.GT.maxits )
369  $ GO TO 100
370 *
371 * Normalize and scale the righthand side vector Pb.
372 *
373  jmax = isamax( blksiz, work( indrv1+1 ), 1 )
374  scl = blksiz*onenrm*max( eps,
375  $ abs( work( indrv4+blksiz ) ) ) /
376  $ abs( work( indrv1+jmax ) )
377  CALL sscal( blksiz, scl, work( indrv1+1 ), 1 )
378 *
379 * Solve the system LU = Pb.
380 *
381  CALL slagts( -1, blksiz, work( indrv4+1 ), work( indrv2+2 ),
382  $ work( indrv3+1 ), work( indrv5+1 ), iwork,
383  $ work( indrv1+1 ), tol, iinfo )
384 *
385 * Reorthogonalize by modified Gram-Schmidt if eigenvalues are
386 * close enough.
387 *
388  IF( jblk.EQ.1 )
389  $ GO TO 90
390  IF( abs( xj-xjm ).GT.ortol )
391  $ gpind = j
392  IF( gpind.NE.j ) THEN
393  DO 80 i = gpind, j - 1
394  ctr = -sdot( blksiz, work( indrv1+1 ), 1, z( b1, i ),
395  $ 1 )
396  CALL saxpy( blksiz, ctr, z( b1, i ), 1,
397  $ work( indrv1+1 ), 1 )
398  80 CONTINUE
399  END IF
400 *
401 * Check the infinity norm of the iterate.
402 *
403  90 CONTINUE
404  jmax = isamax( blksiz, work( indrv1+1 ), 1 )
405  nrm = abs( work( indrv1+jmax ) )
406 *
407 * Continue for additional iterations after norm reaches
408 * stopping criterion.
409 *
410  IF( nrm.LT.stpcrt )
411  $ GO TO 70
412  nrmchk = nrmchk + 1
413  IF( nrmchk.LT.extra+1 )
414  $ GO TO 70
415 *
416  GO TO 110
417 *
418 * If stopping criterion was not satisfied, update info and
419 * store eigenvector number in array ifail.
420 *
421  100 CONTINUE
422  info = info + 1
423  ifail( info ) = j
424 *
425 * Accept iterate as jth eigenvector.
426 *
427  110 CONTINUE
428  scl = one / snrm2( blksiz, work( indrv1+1 ), 1 )
429  jmax = isamax( blksiz, work( indrv1+1 ), 1 )
430  IF( work( indrv1+jmax ).LT.zero )
431  $ scl = -scl
432  CALL sscal( blksiz, scl, work( indrv1+1 ), 1 )
433  120 CONTINUE
434  DO 130 i = 1, n
435  z( i, j ) = zero
436  130 CONTINUE
437  DO 140 i = 1, blksiz
438  z( b1+i-1, j ) = work( indrv1+i )
439  140 CONTINUE
440 *
441 * Save the shift to check eigenvalue spacing at next
442 * iteration.
443 *
444  xjm = xj
445 *
446  150 CONTINUE
447  160 CONTINUE
448 *
449  RETURN
450 *
451 * End of SSTEIN
452 *
real function sdot(N, SX, INCX, SY, INCY)
SDOT
Definition: sdot.f:53
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine slarnv(IDIST, ISEED, N, X)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: slarnv.f:99
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function snrm2(N, X, INCX)
SNRM2
Definition: snrm2.f:56
real function sasum(N, SX, INCX)
SASUM
Definition: sasum.f:54
subroutine saxpy(N, SA, SX, INCX, SY, INCY)
SAXPY
Definition: saxpy.f:54
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine slagtf(N, A, LAMBDA, B, C, TOL, D, IN, INFO)
SLAGTF computes an LU factorization of a matrix T-λI, where T is a general tridiagonal matrix...
Definition: slagtf.f:158
subroutine slagts(JOB, N, A, B, C, D, IN, Y, TOL, INFO)
SLAGTS solves the system of equations (T-λI)x = y or (T-λI)Tx = y,where T is a general tridiagonal ma...
Definition: slagts.f:163
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: