LAPACK  3.10.1
LAPACK: Linear Algebra PACKage

◆ dspevx()

subroutine dspevx ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
double precision, dimension( * )  AP,
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, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

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

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

Purpose:
 DSPEVX 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 DOUBLE PRECISION 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 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 AP 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 "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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, 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 (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 dspevx.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  DOUBLE PRECISION ABSTOL, VL, VU
243 * ..
244 * .. Array Arguments ..
245  INTEGER IFAIL( * ), IWORK( * )
246  DOUBLE PRECISION AP( * ), W( * ), WORK( * ), Z( LDZ, * )
247 * ..
248 *
249 * =====================================================================
250 *
251 * .. Parameters ..
252  DOUBLE PRECISION ZERO, ONE
253  parameter( zero = 0.0d0, one = 1.0d0 )
254 * ..
255 * .. Local Scalars ..
256  LOGICAL ALLEIG, INDEIG, TEST, VALEIG, WANTZ
257  CHARACTER ORDER
258  INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
259  $ INDISP, INDIWO, INDTAU, INDWRK, ISCALE, ITMP1,
260  $ J, JJ, NSPLIT
261  DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
262  $ SIGMA, SMLNUM, TMP1, VLL, VUU
263 * ..
264 * .. External Functions ..
265  LOGICAL LSAME
266  DOUBLE PRECISION DLAMCH, DLANSP
267  EXTERNAL lsame, dlamch, dlansp
268 * ..
269 * .. External Subroutines ..
270  EXTERNAL dcopy, dopgtr, dopmtr, dscal, dsptrd, dstebz,
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( 'DSPEVX', -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 = dlamch( 'Safe minimum' )
341  eps = dlamch( '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  END IF
358  anrm = dlansp( '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 dscal( ( 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 DSPTRD 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 dsptrd( 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 DSTERF or DOPGTR and SSTEQR. If this fails
387 * for some eigenvalue, then try DSTEBZ.
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 dcopy( n, work( indd ), 1, w, 1 )
397  indee = indwrk + 2*n
398  IF( .NOT.wantz ) THEN
399  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
400  CALL dsterf( n, w, work( indee ), info )
401  ELSE
402  CALL dopgtr( uplo, n, ap, work( indtau ), z, ldz,
403  $ work( indwrk ), iinfo )
404  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
405  CALL dsteqr( 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 DSTEBZ and, if eigenvectors are desired, SSTEIN.
421 *
422  IF( wantz ) THEN
423  order = 'B'
424  ELSE
425  order = 'E'
426  END IF
427  indibl = 1
428  indisp = indibl + n
429  indiwo = indisp + n
430  CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
431  $ work( indd ), work( inde ), m, nsplit, w,
432  $ iwork( indibl ), iwork( indisp ), work( indwrk ),
433  $ iwork( indiwo ), info )
434 *
435  IF( wantz ) THEN
436  CALL dstein( n, work( indd ), work( inde ), m, w,
437  $ iwork( indibl ), iwork( indisp ), z, ldz,
438  $ work( indwrk ), iwork( indiwo ), ifail, info )
439 *
440 * Apply orthogonal matrix used in reduction to tridiagonal
441 * form to eigenvectors returned by DSTEIN.
442 *
443  CALL dopmtr( 'L', uplo, 'N', n, m, ap, work( indtau ), z, ldz,
444  $ work( indwrk ), iinfo )
445  END IF
446 *
447 * If matrix was scaled, then rescale eigenvalues appropriately.
448 *
449  20 CONTINUE
450  IF( iscale.EQ.1 ) THEN
451  IF( info.EQ.0 ) THEN
452  imax = m
453  ELSE
454  imax = info - 1
455  END IF
456  CALL dscal( imax, one / sigma, w, 1 )
457  END IF
458 *
459 * If eigenvalues are not in order, then sort them, along with
460 * eigenvectors.
461 *
462  IF( wantz ) THEN
463  DO 40 j = 1, m - 1
464  i = 0
465  tmp1 = w( j )
466  DO 30 jj = j + 1, m
467  IF( w( jj ).LT.tmp1 ) THEN
468  i = jj
469  tmp1 = w( jj )
470  END IF
471  30 CONTINUE
472 *
473  IF( i.NE.0 ) THEN
474  itmp1 = iwork( indibl+i-1 )
475  w( i ) = w( j )
476  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
477  w( j ) = tmp1
478  iwork( indibl+j-1 ) = itmp1
479  CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
480  IF( info.NE.0 ) THEN
481  itmp1 = ifail( i )
482  ifail( i ) = ifail( j )
483  ifail( j ) = itmp1
484  END IF
485  END IF
486  40 CONTINUE
487  END IF
488 *
489  RETURN
490 *
491 * End of DSPEVX
492 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:273
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:131
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:86
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:82
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:79
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:82
double precision function dlansp(NORM, UPLO, N, AP, WORK)
DLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: dlansp.f:114
subroutine dsptrd(UPLO, N, AP, D, E, TAU, INFO)
DSPTRD
Definition: dsptrd.f:150
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
Definition: dstein.f:174
subroutine dopgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
DOPGTR
Definition: dopgtr.f:114
subroutine dopmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
DOPMTR
Definition: dopmtr.f:150
Here is the call graph for this function:
Here is the caller graph for this function: