LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cheevx ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real  VL,
real  VU,
integer  IL,
integer  IU,
real  ABSTOL,
integer  M,
real, dimension( * )  W,
complex, dimension( ldz, * )  Z,
integer  LDZ,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

CHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices

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

Purpose:
 CHEEVX computes selected eigenvalues and, optionally, eigenvectors
 of a complex Hermitian 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 COMPLEX array, dimension (LDA, N)
          On entry, the Hermitian 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 "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)
          On normal exit, the first M elements contain the selected
          eigenvalues in ascending order.
[out]Z
          Z is COMPLEX 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 COMPLEX 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 2*N.
          For optimal efficiency, LWORK >= (NB+1)*N,
          where NB is the max of the blocksize for CHETRD and for
          CUNMTR as 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]RWORK
          RWORK is REAL array, dimension (7*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.
Date
June 2016

Definition at line 261 of file cheevx.f.

261 *
262 * -- LAPACK driver routine (version 3.6.1) --
263 * -- LAPACK is a software package provided by Univ. of Tennessee, --
264 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
265 * June 2016
266 *
267 * .. Scalar Arguments ..
268  CHARACTER jobz, range, uplo
269  INTEGER il, info, iu, lda, ldz, lwork, m, n
270  REAL abstol, vl, vu
271 * ..
272 * .. Array Arguments ..
273  INTEGER ifail( * ), iwork( * )
274  REAL rwork( * ), w( * )
275  COMPLEX a( lda, * ), work( * ), z( ldz, * )
276 * ..
277 *
278 * =====================================================================
279 *
280 * .. Parameters ..
281  REAL zero, one
282  parameter ( zero = 0.0e+0, one = 1.0e+0 )
283  COMPLEX cone
284  parameter ( cone = ( 1.0e+0, 0.0e+0 ) )
285 * ..
286 * .. Local Scalars ..
287  LOGICAL alleig, indeig, lower, lquery, test, valeig,
288  $ wantz
289  CHARACTER order
290  INTEGER i, iinfo, imax, indd, inde, indee, indibl,
291  $ indisp, indiwk, indrwk, indtau, indwrk, iscale,
292  $ itmp1, j, jj, llwork, lwkmin, lwkopt, nb,
293  $ nsplit
294  REAL abstll, anrm, bignum, eps, rmax, rmin, safmin,
295  $ sigma, smlnum, tmp1, vll, vuu
296 * ..
297 * .. External Functions ..
298  LOGICAL lsame
299  INTEGER ilaenv
300  REAL slamch, clanhe
301  EXTERNAL lsame, ilaenv, slamch, clanhe
302 * ..
303 * .. External Subroutines ..
304  EXTERNAL scopy, sscal, sstebz, ssterf, xerbla, csscal,
306  $ cunmtr
307 * ..
308 * .. Intrinsic Functions ..
309  INTRINSIC REAL, max, min, sqrt
310 * ..
311 * .. Executable Statements ..
312 *
313 * Test the input parameters.
314 *
315  lower = lsame( uplo, 'L' )
316  wantz = lsame( jobz, 'V' )
317  alleig = lsame( range, 'A' )
318  valeig = lsame( range, 'V' )
319  indeig = lsame( range, 'I' )
320  lquery = ( lwork.EQ.-1 )
321 *
322  info = 0
323  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
324  info = -1
325  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
326  info = -2
327  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
328  info = -3
329  ELSE IF( n.LT.0 ) THEN
330  info = -4
331  ELSE IF( lda.LT.max( 1, n ) ) THEN
332  info = -6
333  ELSE
334  IF( valeig ) THEN
335  IF( n.GT.0 .AND. vu.LE.vl )
336  $ info = -8
337  ELSE IF( indeig ) THEN
338  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
339  info = -9
340  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
341  info = -10
342  END IF
343  END IF
344  END IF
345  IF( info.EQ.0 ) THEN
346  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
347  info = -15
348  END IF
349  END IF
350 *
351  IF( info.EQ.0 ) THEN
352  IF( n.LE.1 ) THEN
353  lwkmin = 1
354  work( 1 ) = lwkmin
355  ELSE
356  lwkmin = 2*n
357  nb = ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 )
358  nb = max( nb, ilaenv( 1, 'CUNMTR', uplo, n, -1, -1, -1 ) )
359  lwkopt = max( 1, ( nb + 1 )*n )
360  work( 1 ) = lwkopt
361  END IF
362 *
363  IF( lwork.LT.lwkmin .AND. .NOT.lquery )
364  $ info = -17
365  END IF
366 *
367  IF( info.NE.0 ) THEN
368  CALL xerbla( 'CHEEVX', -info )
369  RETURN
370  ELSE IF( lquery ) THEN
371  RETURN
372  END IF
373 *
374 * Quick return if possible
375 *
376  m = 0
377  IF( n.EQ.0 ) THEN
378  RETURN
379  END IF
380 *
381  IF( n.EQ.1 ) THEN
382  IF( alleig .OR. indeig ) THEN
383  m = 1
384  w( 1 ) = a( 1, 1 )
385  ELSE IF( valeig ) THEN
386  IF( vl.LT.REAL( A( 1, 1 ) ) .AND. vu.GE.REAL( A( 1, 1 ) ) )
387  $ THEN
388  m = 1
389  w( 1 ) = a( 1, 1 )
390  END IF
391  END IF
392  IF( wantz )
393  $ z( 1, 1 ) = cone
394  RETURN
395  END IF
396 *
397 * Get machine constants.
398 *
399  safmin = slamch( 'Safe minimum' )
400  eps = slamch( 'Precision' )
401  smlnum = safmin / eps
402  bignum = one / smlnum
403  rmin = sqrt( smlnum )
404  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
405 *
406 * Scale matrix to allowable range, if necessary.
407 *
408  iscale = 0
409  abstll = abstol
410  IF( valeig ) THEN
411  vll = vl
412  vuu = vu
413  END IF
414  anrm = clanhe( 'M', uplo, n, a, lda, rwork )
415  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
416  iscale = 1
417  sigma = rmin / anrm
418  ELSE IF( anrm.GT.rmax ) THEN
419  iscale = 1
420  sigma = rmax / anrm
421  END IF
422  IF( iscale.EQ.1 ) THEN
423  IF( lower ) THEN
424  DO 10 j = 1, n
425  CALL csscal( n-j+1, sigma, a( j, j ), 1 )
426  10 CONTINUE
427  ELSE
428  DO 20 j = 1, n
429  CALL csscal( j, sigma, a( 1, j ), 1 )
430  20 CONTINUE
431  END IF
432  IF( abstol.GT.0 )
433  $ abstll = abstol*sigma
434  IF( valeig ) THEN
435  vll = vl*sigma
436  vuu = vu*sigma
437  END IF
438  END IF
439 *
440 * Call CHETRD to reduce Hermitian matrix to tridiagonal form.
441 *
442  indd = 1
443  inde = indd + n
444  indrwk = inde + n
445  indtau = 1
446  indwrk = indtau + n
447  llwork = lwork - indwrk + 1
448  CALL chetrd( uplo, n, a, lda, rwork( indd ), rwork( inde ),
449  $ work( indtau ), work( indwrk ), llwork, iinfo )
450 *
451 * If all eigenvalues are desired and ABSTOL is less than or equal to
452 * zero, then call SSTERF or CUNGTR and CSTEQR. If this fails for
453 * some eigenvalue, then try SSTEBZ.
454 *
455  test = .false.
456  IF( indeig ) THEN
457  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
458  test = .true.
459  END IF
460  END IF
461  IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
462  CALL scopy( n, rwork( indd ), 1, w, 1 )
463  indee = indrwk + 2*n
464  IF( .NOT.wantz ) THEN
465  CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
466  CALL ssterf( n, w, rwork( indee ), info )
467  ELSE
468  CALL clacpy( 'A', n, n, a, lda, z, ldz )
469  CALL cungtr( uplo, n, z, ldz, work( indtau ),
470  $ work( indwrk ), llwork, iinfo )
471  CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
472  CALL csteqr( jobz, n, w, rwork( indee ), z, ldz,
473  $ rwork( indrwk ), info )
474  IF( info.EQ.0 ) THEN
475  DO 30 i = 1, n
476  ifail( i ) = 0
477  30 CONTINUE
478  END IF
479  END IF
480  IF( info.EQ.0 ) THEN
481  m = n
482  GO TO 40
483  END IF
484  info = 0
485  END IF
486 *
487 * Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
488 *
489  IF( wantz ) THEN
490  order = 'B'
491  ELSE
492  order = 'E'
493  END IF
494  indibl = 1
495  indisp = indibl + n
496  indiwk = indisp + n
497  CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
498  $ rwork( indd ), rwork( inde ), m, nsplit, w,
499  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
500  $ iwork( indiwk ), info )
501 *
502  IF( wantz ) THEN
503  CALL cstein( n, rwork( indd ), rwork( inde ), m, w,
504  $ iwork( indibl ), iwork( indisp ), z, ldz,
505  $ rwork( indrwk ), iwork( indiwk ), ifail, info )
506 *
507 * Apply unitary matrix used in reduction to tridiagonal
508 * form to eigenvectors returned by CSTEIN.
509 *
510  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
511  $ ldz, work( indwrk ), llwork, iinfo )
512  END IF
513 *
514 * If matrix was scaled, then rescale eigenvalues appropriately.
515 *
516  40 CONTINUE
517  IF( iscale.EQ.1 ) THEN
518  IF( info.EQ.0 ) THEN
519  imax = m
520  ELSE
521  imax = info - 1
522  END IF
523  CALL sscal( imax, one / sigma, w, 1 )
524  END IF
525 *
526 * If eigenvalues are not in order, then sort them, along with
527 * eigenvectors.
528 *
529  IF( wantz ) THEN
530  DO 60 j = 1, m - 1
531  i = 0
532  tmp1 = w( j )
533  DO 50 jj = j + 1, m
534  IF( w( jj ).LT.tmp1 ) THEN
535  i = jj
536  tmp1 = w( jj )
537  END IF
538  50 CONTINUE
539 *
540  IF( i.NE.0 ) THEN
541  itmp1 = iwork( indibl+i-1 )
542  w( i ) = w( j )
543  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
544  w( j ) = tmp1
545  iwork( indibl+j-1 ) = itmp1
546  CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
547  IF( info.NE.0 ) THEN
548  itmp1 = ifail( i )
549  ifail( i ) = ifail( j )
550  ifail( j ) = itmp1
551  END IF
552  END IF
553  60 CONTINUE
554  END IF
555 *
556 * Set WORK(1) to optimal complex workspace size.
557 *
558  work( 1 ) = lwkopt
559 *
560  RETURN
561 *
562 * End of CHEEVX
563 *
subroutine sstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
SSTEBZ
Definition: sstebz.f:275
subroutine cunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMTR
Definition: cunmtr.f:174
real function clanhe(NORM, UPLO, N, A, LDA, WORK)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix.
Definition: clanhe.f:126
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine csteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
CSTEQR
Definition: csteqr.f:134
subroutine cstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CSTEIN
Definition: cstein.f:184
subroutine chetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
CHETRD
Definition: chetrd.f:194
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine cungtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
CUNGTR
Definition: cungtr.f:125
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
Definition: cswap.f:52
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: