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

◆ ssyevx()

subroutine ssyevx ( character jobz,
character range,
character uplo,
integer n,
real, dimension( lda, * ) a,
integer lda,
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 lwork,
integer, dimension( * ) iwork,
integer, dimension( * ) ifail,
integer info )

SSYEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices

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

Purpose:
!>
!> SSYEVX computes selected eigenvalues and, optionally, eigenvectors
!> of a real symmetric matrix A.  Eigenvalues and eigenvectors 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]A
!>          A is REAL array, dimension (LDA, N)
!>          On entry, the symmetric matrix A.  If UPLO = 'U', the
!>          leading N-by-N upper triangular part of A contains the
!>          upper triangular part of the matrix A.  If UPLO = 'L',
!>          the leading N-by-N lower triangular part of A contains
!>          the lower triangular part of the matrix A.
!>          On exit, the lower triangle (if UPLO='L') or the upper
!>          triangle (if UPLO='U') of A, including the diagonal, is
!>          destroyed.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[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 A 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  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)
!>          On normal exit, the first M elements contain 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 (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of the array WORK.  LWORK >= 1, when N <= 1;
!>          otherwise 8*N.
!>          For optimal efficiency, LWORK >= (NB+3)*N,
!>          where NB is the max of the blocksize for SSYTRD and SORMTR
!>          returned by ILAENV.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[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 248 of file ssyevx.f.

252*
253* -- LAPACK driver routine --
254* -- LAPACK is a software package provided by Univ. of Tennessee, --
255* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
256*
257* .. Scalar Arguments ..
258 CHARACTER JOBZ, RANGE, UPLO
259 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
260 REAL ABSTOL, VL, VU
261* ..
262* .. Array Arguments ..
263 INTEGER IFAIL( * ), IWORK( * )
264 REAL A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
265* ..
266*
267* =====================================================================
268*
269* .. Parameters ..
270 REAL ZERO, ONE
271 parameter( zero = 0.0e+0, one = 1.0e+0 )
272* ..
273* .. Local Scalars ..
274 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
275 $ WANTZ
276 CHARACTER ORDER
277 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
278 $ INDISP, INDIWO, INDTAU, INDWKN, INDWRK, ISCALE,
279 $ ITMP1, J, JJ, LLWORK, LLWRKN, LWKMIN,
280 $ LWKOPT, NB, NSPLIT
281 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
282 $ SIGMA, SMLNUM, TMP1, VLL, VUU
283* ..
284* .. External Functions ..
285 LOGICAL LSAME
286 INTEGER ILAENV
287 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
288 EXTERNAL lsame, ilaenv, slamch,
290* ..
291* .. External Subroutines ..
292 EXTERNAL scopy, slacpy, sorgtr, sormtr, sscal,
293 $ sstebz,
295* ..
296* .. Intrinsic Functions ..
297 INTRINSIC max, min, sqrt
298* ..
299* .. Executable Statements ..
300*
301* Test the input parameters.
302*
303 lower = lsame( uplo, 'L' )
304 wantz = lsame( jobz, 'V' )
305 alleig = lsame( range, 'A' )
306 valeig = lsame( range, 'V' )
307 indeig = lsame( range, 'I' )
308 lquery = ( lwork.EQ.-1 )
309*
310 info = 0
311 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
312 info = -1
313 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
314 info = -2
315 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
316 info = -3
317 ELSE IF( n.LT.0 ) THEN
318 info = -4
319 ELSE IF( lda.LT.max( 1, n ) ) THEN
320 info = -6
321 ELSE
322 IF( valeig ) THEN
323 IF( n.GT.0 .AND. vu.LE.vl )
324 $ info = -8
325 ELSE IF( indeig ) THEN
326 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
327 info = -9
328 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
329 info = -10
330 END IF
331 END IF
332 END IF
333 IF( info.EQ.0 ) THEN
334 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
335 info = -15
336 END IF
337 END IF
338*
339 IF( info.EQ.0 ) THEN
340 IF( n.LE.1 ) THEN
341 lwkmin = 1
342 lwkopt = 1
343 ELSE
344 lwkmin = 8*n
345 nb = ilaenv( 1, 'SSYTRD', uplo, n, -1, -1, -1 )
346 nb = max( nb, ilaenv( 1, 'SORMTR', uplo, n, -1, -1,
347 $ -1 ) )
348 lwkopt = max( lwkmin, ( nb + 3 )*n )
349 END IF
350 work( 1 ) = sroundup_lwork( lwkopt )
351*
352 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
353 $ info = -17
354 END IF
355*
356 IF( info.NE.0 ) THEN
357 CALL xerbla( 'SSYEVX', -info )
358 RETURN
359 ELSE IF( lquery ) THEN
360 RETURN
361 END IF
362*
363* Quick return if possible
364*
365 m = 0
366 IF( n.EQ.0 ) THEN
367 RETURN
368 END IF
369*
370 IF( n.EQ.1 ) THEN
371 IF( alleig .OR. indeig ) THEN
372 m = 1
373 w( 1 ) = a( 1, 1 )
374 ELSE
375 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
376 m = 1
377 w( 1 ) = a( 1, 1 )
378 END IF
379 END IF
380 IF( wantz )
381 $ z( 1, 1 ) = one
382 RETURN
383 END IF
384*
385* Get machine constants.
386*
387 safmin = slamch( 'Safe minimum' )
388 eps = slamch( 'Precision' )
389 smlnum = safmin / eps
390 bignum = one / smlnum
391 rmin = sqrt( smlnum )
392 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
393*
394* Scale matrix to allowable range, if necessary.
395*
396 iscale = 0
397 abstll = abstol
398 IF( valeig ) THEN
399 vll = vl
400 vuu = vu
401 END IF
402 anrm = slansy( 'M', uplo, n, a, lda, work )
403 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
404 iscale = 1
405 sigma = rmin / anrm
406 ELSE IF( anrm.GT.rmax ) THEN
407 iscale = 1
408 sigma = rmax / anrm
409 END IF
410 IF( iscale.EQ.1 ) THEN
411 IF( lower ) THEN
412 DO 10 j = 1, n
413 CALL sscal( n-j+1, sigma, a( j, j ), 1 )
414 10 CONTINUE
415 ELSE
416 DO 20 j = 1, n
417 CALL sscal( j, sigma, a( 1, j ), 1 )
418 20 CONTINUE
419 END IF
420 IF( abstol.GT.0 )
421 $ abstll = abstol*sigma
422 IF( valeig ) THEN
423 vll = vl*sigma
424 vuu = vu*sigma
425 END IF
426 END IF
427*
428* Call SSYTRD to reduce symmetric matrix to tridiagonal form.
429*
430 indtau = 1
431 inde = indtau + n
432 indd = inde + n
433 indwrk = indd + n
434 llwork = lwork - indwrk + 1
435 CALL ssytrd( uplo, n, a, lda, work( indd ), work( inde ),
436 $ work( indtau ), work( indwrk ), llwork, iinfo )
437*
438* If all eigenvalues are desired and ABSTOL is less than or equal to
439* zero, then call SSTERF or SORGTR and SSTEQR. If this fails for
440* some eigenvalue, then try SSTEBZ.
441*
442 test = .false.
443 IF( indeig ) THEN
444 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
445 test = .true.
446 END IF
447 END IF
448 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
449 CALL scopy( n, work( indd ), 1, w, 1 )
450 indee = indwrk + 2*n
451 IF( .NOT.wantz ) THEN
452 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
453 CALL ssterf( n, w, work( indee ), info )
454 ELSE
455 CALL slacpy( 'A', n, n, a, lda, z, ldz )
456 CALL sorgtr( uplo, n, z, ldz, work( indtau ),
457 $ work( indwrk ), llwork, iinfo )
458 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
459 CALL ssteqr( jobz, n, w, work( indee ), z, ldz,
460 $ work( indwrk ), info )
461 IF( info.EQ.0 ) THEN
462 DO 30 i = 1, n
463 ifail( i ) = 0
464 30 CONTINUE
465 END IF
466 END IF
467 IF( info.EQ.0 ) THEN
468 m = n
469 GO TO 40
470 END IF
471 info = 0
472 END IF
473*
474* Otherwise, call SSTEBZ and, if eigenvectors are desired, SSTEIN.
475*
476 IF( wantz ) THEN
477 order = 'B'
478 ELSE
479 order = 'E'
480 END IF
481 indibl = 1
482 indisp = indibl + n
483 indiwo = indisp + n
484 CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
485 $ work( indd ), work( inde ), m, nsplit, w,
486 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
487 $ iwork( indiwo ), info )
488*
489 IF( wantz ) THEN
490 CALL sstein( n, work( indd ), work( inde ), m, w,
491 $ iwork( indibl ), iwork( indisp ), z, ldz,
492 $ work( indwrk ), iwork( indiwo ), ifail, info )
493*
494* Apply orthogonal matrix used in reduction to tridiagonal
495* form to eigenvectors returned by SSTEIN.
496*
497 indwkn = inde
498 llwrkn = lwork - indwkn + 1
499 CALL sormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
500 $ z,
501 $ ldz, work( indwkn ), llwrkn, iinfo )
502 END IF
503*
504* If matrix was scaled, then rescale eigenvalues appropriately.
505*
506 40 CONTINUE
507 IF( iscale.EQ.1 ) THEN
508 IF( info.EQ.0 ) THEN
509 imax = m
510 ELSE
511 imax = info - 1
512 END IF
513 CALL sscal( imax, one / sigma, w, 1 )
514 END IF
515*
516* If eigenvalues are not in order, then sort them, along with
517* eigenvectors.
518*
519 IF( wantz ) THEN
520 DO 60 j = 1, m - 1
521 i = 0
522 tmp1 = w( j )
523 DO 50 jj = j + 1, m
524 IF( w( jj ).LT.tmp1 ) THEN
525 i = jj
526 tmp1 = w( jj )
527 END IF
528 50 CONTINUE
529*
530 IF( i.NE.0 ) THEN
531 itmp1 = iwork( indibl+i-1 )
532 w( i ) = w( j )
533 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
534 w( j ) = tmp1
535 iwork( indibl+j-1 ) = itmp1
536 CALL sswap( n, z( 1, i ), 1, z( 1, j ), 1 )
537 IF( info.NE.0 ) THEN
538 itmp1 = ifail( i )
539 ifail( i ) = ifail( j )
540 ifail( j ) = itmp1
541 END IF
542 END IF
543 60 CONTINUE
544 END IF
545*
546* Set WORK(1) to optimal workspace size.
547*
548 work( 1 ) = sroundup_lwork( lwkopt )
549*
550 RETURN
551*
552* End of SSYEVX
553*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine ssytrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
SSYTRD
Definition ssytrd.f:191
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:120
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
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:272
subroutine sstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
SSTEIN
Definition sstein.f:172
subroutine ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:129
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine sorgtr(uplo, n, a, lda, tau, work, lwork, info)
SORGTR
Definition sorgtr.f:121
subroutine sormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
SORMTR
Definition sormtr.f:171
Here is the call graph for this function:
Here is the caller graph for this function: