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

◆ dsyevx()

subroutine dsyevx ( character jobz,
character range,
character uplo,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
double precision vl,
double precision vu,
integer il,
integer iu,
double precision abstol,
integer m,
double precision, dimension( * ) w,
double precision, dimension( ldz, * ) z,
integer ldz,
double precision, dimension( * ) work,
integer lwork,
integer, dimension( * ) iwork,
integer, dimension( * ) ifail,
integer info )

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

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

Purpose:
!>
!> DSYEVX 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 DOUBLE PRECISION 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 DOUBLE PRECISION
!>          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 DOUBLE PRECISION
!>          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 DOUBLE PRECISION
!>          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*DLAMCH('S'), not zero.
!>          If this routine returns with INFO>0, indicating that some
!>          eigenvectors did not converge, try setting ABSTOL to
!>          2*DLAMCH('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 DOUBLE PRECISION array, dimension (N)
!>          On normal exit, the first M elements contain the selected
!>          eigenvalues in ascending order.
!> 
[out]Z
!>          Z is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DSYTRD and DORMTR
!>          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 dsyevx.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 DOUBLE PRECISION ABSTOL, VL, VU
261* ..
262* .. Array Arguments ..
263 INTEGER IFAIL( * ), IWORK( * )
264 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
265* ..
266*
267* =====================================================================
268*
269* .. Parameters ..
270 DOUBLE PRECISION ZERO, ONE
271 parameter( zero = 0.0d+0, one = 1.0d+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 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
282 $ SIGMA, SMLNUM, TMP1, VLL, VUU
283* ..
284* .. External Functions ..
285 LOGICAL LSAME
286 INTEGER ILAENV
287 DOUBLE PRECISION DLAMCH, DLANSY
288 EXTERNAL lsame, ilaenv, dlamch, dlansy
289* ..
290* .. External Subroutines ..
291 EXTERNAL dcopy, dlacpy, dorgtr, dormtr, dscal,
292 $ dstebz,
294* ..
295* .. Intrinsic Functions ..
296 INTRINSIC max, min, sqrt
297* ..
298* .. Executable Statements ..
299*
300* Test the input parameters.
301*
302 lower = lsame( uplo, 'L' )
303 wantz = lsame( jobz, 'V' )
304 alleig = lsame( range, 'A' )
305 valeig = lsame( range, 'V' )
306 indeig = lsame( range, 'I' )
307 lquery = ( lwork.EQ.-1 )
308*
309 info = 0
310 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
311 info = -1
312 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
313 info = -2
314 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
315 info = -3
316 ELSE IF( n.LT.0 ) THEN
317 info = -4
318 ELSE IF( lda.LT.max( 1, n ) ) THEN
319 info = -6
320 ELSE
321 IF( valeig ) THEN
322 IF( n.GT.0 .AND. vu.LE.vl )
323 $ info = -8
324 ELSE IF( indeig ) THEN
325 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
326 info = -9
327 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
328 info = -10
329 END IF
330 END IF
331 END IF
332 IF( info.EQ.0 ) THEN
333 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
334 info = -15
335 END IF
336 END IF
337*
338 IF( info.EQ.0 ) THEN
339 IF( n.LE.1 ) THEN
340 lwkmin = 1
341 lwkopt = 1
342 ELSE
343 lwkmin = 8*n
344 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
345 nb = max( nb, ilaenv( 1, 'DORMTR', uplo, n, -1, -1,
346 $ -1 ) )
347 lwkopt = max( lwkmin, ( nb + 3 )*n )
348 END IF
349 work( 1 ) = lwkopt
350*
351 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
352 $ info = -17
353 END IF
354*
355 IF( info.NE.0 ) THEN
356 CALL xerbla( 'DSYEVX', -info )
357 RETURN
358 ELSE IF( lquery ) THEN
359 RETURN
360 END IF
361*
362* Quick return if possible
363*
364 m = 0
365 IF( n.EQ.0 ) THEN
366 RETURN
367 END IF
368*
369 IF( n.EQ.1 ) THEN
370 IF( alleig .OR. indeig ) THEN
371 m = 1
372 w( 1 ) = a( 1, 1 )
373 ELSE
374 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
375 m = 1
376 w( 1 ) = a( 1, 1 )
377 END IF
378 END IF
379 IF( wantz )
380 $ z( 1, 1 ) = one
381 RETURN
382 END IF
383*
384* Get machine constants.
385*
386 safmin = dlamch( 'Safe minimum' )
387 eps = dlamch( 'Precision' )
388 smlnum = safmin / eps
389 bignum = one / smlnum
390 rmin = sqrt( smlnum )
391 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
392*
393* Scale matrix to allowable range, if necessary.
394*
395 iscale = 0
396 abstll = abstol
397 IF( valeig ) THEN
398 vll = vl
399 vuu = vu
400 END IF
401 anrm = dlansy( 'M', uplo, n, a, lda, work )
402 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
403 iscale = 1
404 sigma = rmin / anrm
405 ELSE IF( anrm.GT.rmax ) THEN
406 iscale = 1
407 sigma = rmax / anrm
408 END IF
409 IF( iscale.EQ.1 ) THEN
410 IF( lower ) THEN
411 DO 10 j = 1, n
412 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
413 10 CONTINUE
414 ELSE
415 DO 20 j = 1, n
416 CALL dscal( j, sigma, a( 1, j ), 1 )
417 20 CONTINUE
418 END IF
419 IF( abstol.GT.0 )
420 $ abstll = abstol*sigma
421 IF( valeig ) THEN
422 vll = vl*sigma
423 vuu = vu*sigma
424 END IF
425 END IF
426*
427* Call DSYTRD to reduce symmetric matrix to tridiagonal form.
428*
429 indtau = 1
430 inde = indtau + n
431 indd = inde + n
432 indwrk = indd + n
433 llwork = lwork - indwrk + 1
434 CALL dsytrd( uplo, n, a, lda, work( indd ), work( inde ),
435 $ work( indtau ), work( indwrk ), llwork, iinfo )
436*
437* If all eigenvalues are desired and ABSTOL is less than or equal to
438* zero, then call DSTERF or DORGTR and SSTEQR. If this fails for
439* some eigenvalue, then try DSTEBZ.
440*
441 test = .false.
442 IF( indeig ) THEN
443 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
444 test = .true.
445 END IF
446 END IF
447 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
448 CALL dcopy( n, work( indd ), 1, w, 1 )
449 indee = indwrk + 2*n
450 IF( .NOT.wantz ) THEN
451 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
452 CALL dsterf( n, w, work( indee ), info )
453 ELSE
454 CALL dlacpy( 'A', n, n, a, lda, z, ldz )
455 CALL dorgtr( uplo, n, z, ldz, work( indtau ),
456 $ work( indwrk ), llwork, iinfo )
457 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
458 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
459 $ work( indwrk ), info )
460 IF( info.EQ.0 ) THEN
461 DO 30 i = 1, n
462 ifail( i ) = 0
463 30 CONTINUE
464 END IF
465 END IF
466 IF( info.EQ.0 ) THEN
467 m = n
468 GO TO 40
469 END IF
470 info = 0
471 END IF
472*
473* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
474*
475 IF( wantz ) THEN
476 order = 'B'
477 ELSE
478 order = 'E'
479 END IF
480 indibl = 1
481 indisp = indibl + n
482 indiwo = indisp + n
483 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
484 $ work( indd ), work( inde ), m, nsplit, w,
485 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
486 $ iwork( indiwo ), info )
487*
488 IF( wantz ) THEN
489 CALL dstein( n, work( indd ), work( inde ), m, w,
490 $ iwork( indibl ), iwork( indisp ), z, ldz,
491 $ work( indwrk ), iwork( indiwo ), ifail, info )
492*
493* Apply orthogonal matrix used in reduction to tridiagonal
494* form to eigenvectors returned by DSTEIN.
495*
496 indwkn = inde
497 llwrkn = lwork - indwkn + 1
498 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
499 $ z,
500 $ ldz, work( indwkn ), llwrkn, iinfo )
501 END IF
502*
503* If matrix was scaled, then rescale eigenvalues appropriately.
504*
505 40 CONTINUE
506 IF( iscale.EQ.1 ) THEN
507 IF( info.EQ.0 ) THEN
508 imax = m
509 ELSE
510 imax = info - 1
511 END IF
512 CALL dscal( imax, one / sigma, w, 1 )
513 END IF
514*
515* If eigenvalues are not in order, then sort them, along with
516* eigenvectors.
517*
518 IF( wantz ) THEN
519 DO 60 j = 1, m - 1
520 i = 0
521 tmp1 = w( j )
522 DO 50 jj = j + 1, m
523 IF( w( jj ).LT.tmp1 ) THEN
524 i = jj
525 tmp1 = w( jj )
526 END IF
527 50 CONTINUE
528*
529 IF( i.NE.0 ) THEN
530 itmp1 = iwork( indibl+i-1 )
531 w( i ) = w( j )
532 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
533 w( j ) = tmp1
534 iwork( indibl+j-1 ) = itmp1
535 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
536 IF( info.NE.0 ) THEN
537 itmp1 = ifail( i )
538 ifail( i ) = ifail( j )
539 ifail( j ) = itmp1
540 END IF
541 END IF
542 60 CONTINUE
543 END IF
544*
545* Set WORK(1) to optimal workspace size.
546*
547 work( 1 ) = lwkopt
548*
549 RETURN
550*
551* End of DSYEVX
552*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsytrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
DSYTRD
Definition dsytrd.f:191
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlansy(norm, uplo, n, a, lda, work)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansy.f:120
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:272
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:172
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:129
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dorgtr(uplo, n, a, lda, tau, work, lwork, info)
DORGTR
Definition dorgtr.f:121
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:170
Here is the call graph for this function:
Here is the caller graph for this function: