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

◆ sspevx()

subroutine sspevx ( character  jobz,
character  range,
character  uplo,
integer  n,
real, dimension( * )  ap,
real  vl,
real  vu,
integer  il,
integer  iu,
real  abstol,
integer  m,
real, dimension( * )  w,
real, dimension( ldz, * )  z,
integer  ldz,
real, dimension( * )  work,
integer, dimension( * )  iwork,
integer, dimension( * )  ifail,
integer  info 
)

SSPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

Purpose:
 SSPEVX computes selected eigenvalues and, optionally, eigenvectors
 of a real symmetric matrix A in packed storage.  Eigenvalues/vectors
 can be selected by specifying either a range of values or a range of
 indices for the desired eigenvalues.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]RANGE
          RANGE is CHARACTER*1
          = 'A': all eigenvalues will be found;
          = 'V': all eigenvalues in the half-open interval (VL,VU]
                 will be found;
          = 'I': the IL-th through IU-th eigenvalues will be found.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]AP
          AP is REAL array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, AP is overwritten by values generated during the
          reduction to tridiagonal form.  If UPLO = 'U', the diagonal
          and first superdiagonal of the tridiagonal matrix T overwrite
          the corresponding elements of A, and if UPLO = 'L', the
          diagonal and first subdiagonal of T overwrite the
          corresponding elements of A.
[in]VL
          VL is REAL
          If RANGE='V', the lower bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]VU
          VU is REAL
          If RANGE='V', the upper bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]IL
          IL is INTEGER
          If RANGE='I', the index of the
          smallest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]IU
          IU is INTEGER
          If RANGE='I', the index of the
          largest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]ABSTOL
          ABSTOL is REAL
          The absolute error tolerance for the eigenvalues.
          An approximate eigenvalue is accepted as converged
          when it is determined to lie in an interval [a,b]
          of width less than or equal to

                  ABSTOL + EPS *   max( |a|,|b| ) ,

          where EPS is the machine precision.  If ABSTOL is less than
          or equal to zero, then  EPS*|T|  will be used in its place,
          where |T| is the 1-norm of the tridiagonal matrix obtained
          by reducing AP to tridiagonal form.

          Eigenvalues will be computed most accurately when ABSTOL is
          set to twice the underflow threshold 2*SLAMCH('S'), not zero.
          If this routine returns with INFO>0, indicating that some
          eigenvectors did not converge, try setting ABSTOL to
          2*SLAMCH('S').

          See "Computing Small Singular Values of Bidiagonal Matrices
          with Guaranteed High Relative Accuracy," by Demmel and
          Kahan, LAPACK Working Note #3.
[out]M
          M is INTEGER
          The total number of eigenvalues found.  0 <= M <= N.
          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the selected eigenvalues in ascending order.
[out]Z
          Z is REAL array, dimension (LDZ, max(1,M))
          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
          contain the orthonormal eigenvectors of the matrix A
          corresponding to the selected eigenvalues, with the i-th
          column of Z holding the eigenvector associated with W(i).
          If an eigenvector fails to converge, then that column of Z
          contains the latest approximation to the eigenvector, and the
          index of the eigenvector is returned in IFAIL.
          If JOBZ = 'N', then Z is not referenced.
          Note: the user must ensure that at least max(1,M) columns are
          supplied in the array Z; if RANGE = 'V', the exact value of M
          is not known in advance and an upper bound must be used.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is REAL array, dimension (8*N)
[out]IWORK
          IWORK is INTEGER array, dimension (5*N)
[out]IFAIL
          IFAIL is INTEGER array, dimension (N)
          If JOBZ = 'V', then if INFO = 0, the first M elements of
          IFAIL are zero.  If INFO > 0, then IFAIL contains the
          indices of the eigenvectors that failed to converge.
          If JOBZ = 'N', then IFAIL is not referenced.
[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.
                Their indices are stored in array IFAIL.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 231 of file sspevx.f.

234*
235* -- LAPACK driver routine --
236* -- LAPACK is a software package provided by Univ. of Tennessee, --
237* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
238*
239* .. Scalar Arguments ..
240 CHARACTER JOBZ, RANGE, UPLO
241 INTEGER IL, INFO, IU, LDZ, M, N
242 REAL ABSTOL, VL, VU
243* ..
244* .. Array Arguments ..
245 INTEGER IFAIL( * ), IWORK( * )
246 REAL AP( * ), W( * ), WORK( * ), Z( LDZ, * )
247* ..
248*
249* =====================================================================
250*
251* .. Parameters ..
252 REAL ZERO, ONE
253 parameter( zero = 0.0e0, one = 1.0e0 )
254* ..
255* .. Local Scalars ..
256 LOGICAL ALLEIG, INDEIG, TEST, VALEIG, WANTZ
257 CHARACTER ORDER
258 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE,
259 $ INDISP, INDIWO, INDTAU, INDWRK, ISCALE, ITMP1,
260 $ J, JJ, NSPLIT
261 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
262 $ SIGMA, SMLNUM, TMP1, VLL, VUU
263* ..
264* .. External Functions ..
265 LOGICAL LSAME
266 REAL SLAMCH, SLANSP
267 EXTERNAL lsame, slamch, slansp
268* ..
269* .. External Subroutines ..
270 EXTERNAL scopy, sopgtr, sopmtr, sscal, ssptrd, sstebz,
272* ..
273* .. Intrinsic Functions ..
274 INTRINSIC max, min, sqrt
275* ..
276* .. Executable Statements ..
277*
278* Test the input parameters.
279*
280 wantz = lsame( jobz, 'V' )
281 alleig = lsame( range, 'A' )
282 valeig = lsame( range, 'V' )
283 indeig = lsame( range, 'I' )
284*
285 info = 0
286 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
287 info = -1
288 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
289 info = -2
290 ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
291 $ THEN
292 info = -3
293 ELSE IF( n.LT.0 ) THEN
294 info = -4
295 ELSE
296 IF( valeig ) THEN
297 IF( n.GT.0 .AND. vu.LE.vl )
298 $ info = -7
299 ELSE IF( indeig ) THEN
300 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
301 info = -8
302 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
303 info = -9
304 END IF
305 END IF
306 END IF
307 IF( info.EQ.0 ) THEN
308 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
309 $ info = -14
310 END IF
311*
312 IF( info.NE.0 ) THEN
313 CALL xerbla( 'SSPEVX', -info )
314 RETURN
315 END IF
316*
317* Quick return if possible
318*
319 m = 0
320 IF( n.EQ.0 )
321 $ RETURN
322*
323 IF( n.EQ.1 ) THEN
324 IF( alleig .OR. indeig ) THEN
325 m = 1
326 w( 1 ) = ap( 1 )
327 ELSE
328 IF( vl.LT.ap( 1 ) .AND. vu.GE.ap( 1 ) ) THEN
329 m = 1
330 w( 1 ) = ap( 1 )
331 END IF
332 END IF
333 IF( wantz )
334 $ z( 1, 1 ) = one
335 RETURN
336 END IF
337*
338* Get machine constants.
339*
340 safmin = slamch( 'Safe minimum' )
341 eps = slamch( 'Precision' )
342 smlnum = safmin / eps
343 bignum = one / smlnum
344 rmin = sqrt( smlnum )
345 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
346*
347* Scale matrix to allowable range, if necessary.
348*
349 iscale = 0
350 abstll = abstol
351 IF ( valeig ) THEN
352 vll = vl
353 vuu = vu
354 ELSE
355 vll = zero
356 vuu = zero
357 ENDIF
358 anrm = slansp( 'M', uplo, n, ap, work )
359 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
360 iscale = 1
361 sigma = rmin / anrm
362 ELSE IF( anrm.GT.rmax ) THEN
363 iscale = 1
364 sigma = rmax / anrm
365 END IF
366 IF( iscale.EQ.1 ) THEN
367 CALL sscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
368 IF( abstol.GT.0 )
369 $ abstll = abstol*sigma
370 IF( valeig ) THEN
371 vll = vl*sigma
372 vuu = vu*sigma
373 END IF
374 END IF
375*
376* Call SSPTRD to reduce symmetric packed matrix to tridiagonal form.
377*
378 indtau = 1
379 inde = indtau + n
380 indd = inde + n
381 indwrk = indd + n
382 CALL ssptrd( uplo, n, ap, work( indd ), work( inde ),
383 $ work( indtau ), iinfo )
384*
385* If all eigenvalues are desired and ABSTOL is less than or equal
386* to zero, then call SSTERF or SOPGTR and SSTEQR. If this fails
387* for some eigenvalue, then try SSTEBZ.
388*
389 test = .false.
390 IF (indeig) THEN
391 IF (il.EQ.1 .AND. iu.EQ.n) THEN
392 test = .true.
393 END IF
394 END IF
395 IF ((alleig .OR. test) .AND. (abstol.LE.zero)) THEN
396 CALL scopy( n, work( indd ), 1, w, 1 )
397 indee = indwrk + 2*n
398 IF( .NOT.wantz ) THEN
399 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
400 CALL ssterf( n, w, work( indee ), info )
401 ELSE
402 CALL sopgtr( uplo, n, ap, work( indtau ), z, ldz,
403 $ work( indwrk ), iinfo )
404 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
405 CALL ssteqr( jobz, n, w, work( indee ), z, ldz,
406 $ work( indwrk ), info )
407 IF( info.EQ.0 ) THEN
408 DO 10 i = 1, n
409 ifail( i ) = 0
410 10 CONTINUE
411 END IF
412 END IF
413 IF( info.EQ.0 ) THEN
414 m = n
415 GO TO 20
416 END IF
417 info = 0
418 END IF
419*
420* Otherwise, call SSTEBZ and, if eigenvectors are desired, SSTEIN.
421*
422 IF( wantz ) THEN
423 order = 'B'
424 ELSE
425 order = 'E'
426 END IF
427 indisp = 1 + n
428 indiwo = indisp + n
429 CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
430 $ work( indd ), work( inde ), m, nsplit, w,
431 $ iwork( 1 ), iwork( indisp ), work( indwrk ),
432 $ iwork( indiwo ), info )
433*
434 IF( wantz ) THEN
435 CALL sstein( n, work( indd ), work( inde ), m, w,
436 $ iwork( 1 ), iwork( indisp ), z, ldz,
437 $ work( indwrk ), iwork( indiwo ), ifail, info )
438*
439* Apply orthogonal matrix used in reduction to tridiagonal
440* form to eigenvectors returned by SSTEIN.
441*
442 CALL sopmtr( 'L', uplo, 'N', n, m, ap, work( indtau ), z, ldz,
443 $ work( indwrk ), iinfo )
444 END IF
445*
446* If matrix was scaled, then rescale eigenvalues appropriately.
447*
448 20 CONTINUE
449 IF( iscale.EQ.1 ) THEN
450 IF( info.EQ.0 ) THEN
451 imax = m
452 ELSE
453 imax = info - 1
454 END IF
455 CALL sscal( imax, one / sigma, w, 1 )
456 END IF
457*
458* If eigenvalues are not in order, then sort them, along with
459* eigenvectors.
460*
461 IF( wantz ) THEN
462 DO 40 j = 1, m - 1
463 i = 0
464 tmp1 = w( j )
465 DO 30 jj = j + 1, m
466 IF( w( jj ).LT.tmp1 ) THEN
467 i = jj
468 tmp1 = w( jj )
469 END IF
470 30 CONTINUE
471*
472 IF( i.NE.0 ) THEN
473 itmp1 = iwork( 1 + i-1 )
474 w( i ) = w( j )
475 iwork( 1 + i-1 ) = iwork( 1 + j-1 )
476 w( j ) = tmp1
477 iwork( 1 + j-1 ) = itmp1
478 CALL sswap( n, z( 1, i ), 1, z( 1, j ), 1 )
479 IF( info.NE.0 ) THEN
480 itmp1 = ifail( i )
481 ifail( i ) = ifail( j )
482 ifail( j ) = itmp1
483 END IF
484 END IF
485 40 CONTINUE
486 END IF
487*
488 RETURN
489*
490* End of SSPEVX
491*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine ssptrd(uplo, n, ap, d, e, tau, info)
SSPTRD
Definition ssptrd.f:150
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slansp(norm, uplo, n, ap, work)
SLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansp.f:114
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
SSTEBZ
Definition sstebz.f:273
subroutine sstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
SSTEIN
Definition sstein.f:174
subroutine ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:131
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine sopgtr(uplo, n, ap, tau, q, ldq, work, info)
SOPGTR
Definition sopgtr.f:114
subroutine sopmtr(side, uplo, trans, m, n, ap, tau, c, ldc, work, info)
SOPMTR
Definition sopmtr.f:150
Here is the call graph for this function:
Here is the caller graph for this function: