LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ cheevx_2stage()

subroutine cheevx_2stage ( 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_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices

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

Purpose:
 CHEEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
 of a complex Hermitian matrix A using the 2stage technique for
 the reduction to tridiagonal.  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.
                  Not available in this release.
[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  
          If JOBZ = 'N' and N > 1, LWORK must be queried.
                                   LWORK = MAX(1, 8*N, dimension) where
                                   dimension = max(stage1,stage2) + (KD+1)*N + N
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + N
                                   where KD is the blocking size of the reduction,
                                   FACTOPTNB is the blocking used by the QR or LQ
                                   algorithm, usually FACTOPTNB=128 is a good choice
                                   NTHREADS is the number of threads used when
                                   openMP compilation is enabled, otherwise =1.
          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available

          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.
Further Details:
  All details about the 2stage techniques are available in:

  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
  Parallel reduction to condensed forms for symmetric eigenvalue problems
  using aggregated fine-grained and memory-aware kernels. In Proceedings
  of 2011 International Conference for High Performance Computing,
  Networking, Storage and Analysis (SC '11), New York, NY, USA,
  Article 8 , 11 pages.
  http://doi.acm.org/10.1145/2063384.2063394

  A. Haidar, J. Kurzak, P. Luszczek, 2013.
  An improved parallel singular value algorithm and its implementation 
  for multicore hardware, In Proceedings of 2013 International Conference
  for High Performance Computing, Networking, Storage and Analysis (SC '13).
  Denver, Colorado, USA, 2013.
  Article 90, 12 pages.
  http://doi.acm.org/10.1145/2503210.2503292

  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
  A novel hybrid CPU-GPU generalized eigensolver for electronic structure 
  calculations based on fine-grained memory aware tasks.
  International Journal of High Performance Computing Applications.
  Volume 28 Issue 2, Pages 196-209, May 2014.
  http://hpc.sagepub.com/content/28/2/196 

Definition at line 303 of file cheevx_2stage.f.

306 *
307  IMPLICIT NONE
308 *
309 * -- LAPACK driver routine --
310 * -- LAPACK is a software package provided by Univ. of Tennessee, --
311 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
312 *
313 * .. Scalar Arguments ..
314  CHARACTER JOBZ, RANGE, UPLO
315  INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
316  REAL ABSTOL, VL, VU
317 * ..
318 * .. Array Arguments ..
319  INTEGER IFAIL( * ), IWORK( * )
320  REAL RWORK( * ), W( * )
321  COMPLEX A( LDA, * ), WORK( * ), Z( LDZ, * )
322 * ..
323 *
324 * =====================================================================
325 *
326 * .. Parameters ..
327  REAL ZERO, ONE
328  parameter( zero = 0.0e+0, one = 1.0e+0 )
329  COMPLEX CONE
330  parameter( cone = ( 1.0e+0, 0.0e+0 ) )
331 * ..
332 * .. Local Scalars ..
333  LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
334  $ WANTZ
335  CHARACTER ORDER
336  INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
337  $ INDISP, INDIWK, INDRWK, INDTAU, INDWRK, ISCALE,
338  $ ITMP1, J, JJ, LLWORK,
339  $ NSPLIT, LWMIN, LHTRD, LWTRD, KD, IB, INDHOUS
340  REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
341  $ SIGMA, SMLNUM, TMP1, VLL, VUU
342 * ..
343 * .. External Functions ..
344  LOGICAL LSAME
345  INTEGER ILAENV2STAGE
346  REAL SLAMCH, CLANHE
347  EXTERNAL lsame, slamch, clanhe, ilaenv2stage
348 * ..
349 * .. External Subroutines ..
350  EXTERNAL scopy, sscal, sstebz, ssterf, xerbla, csscal,
352  $ chetrd_2stage
353 * ..
354 * .. Intrinsic Functions ..
355  INTRINSIC real, max, min, sqrt
356 * ..
357 * .. Executable Statements ..
358 *
359 * Test the input parameters.
360 *
361  lower = lsame( uplo, 'L' )
362  wantz = lsame( jobz, 'V' )
363  alleig = lsame( range, 'A' )
364  valeig = lsame( range, 'V' )
365  indeig = lsame( range, 'I' )
366  lquery = ( lwork.EQ.-1 )
367 *
368  info = 0
369  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
370  info = -1
371  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
372  info = -2
373  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
374  info = -3
375  ELSE IF( n.LT.0 ) THEN
376  info = -4
377  ELSE IF( lda.LT.max( 1, n ) ) THEN
378  info = -6
379  ELSE
380  IF( valeig ) THEN
381  IF( n.GT.0 .AND. vu.LE.vl )
382  $ info = -8
383  ELSE IF( indeig ) THEN
384  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
385  info = -9
386  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
387  info = -10
388  END IF
389  END IF
390  END IF
391  IF( info.EQ.0 ) THEN
392  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
393  info = -15
394  END IF
395  END IF
396 *
397  IF( info.EQ.0 ) THEN
398  IF( n.LE.1 ) THEN
399  lwmin = 1
400  work( 1 ) = lwmin
401  ELSE
402  kd = ilaenv2stage( 1, 'CHETRD_2STAGE', jobz,
403  $ n, -1, -1, -1 )
404  ib = ilaenv2stage( 2, 'CHETRD_2STAGE', jobz,
405  $ n, kd, -1, -1 )
406  lhtrd = ilaenv2stage( 3, 'CHETRD_2STAGE', jobz,
407  $ n, kd, ib, -1 )
408  lwtrd = ilaenv2stage( 4, 'CHETRD_2STAGE', jobz,
409  $ n, kd, ib, -1 )
410  lwmin = n + lhtrd + lwtrd
411  work( 1 ) = lwmin
412  END IF
413 *
414  IF( lwork.LT.lwmin .AND. .NOT.lquery )
415  $ info = -17
416  END IF
417 *
418  IF( info.NE.0 ) THEN
419  CALL xerbla( 'CHEEVX_2STAGE', -info )
420  RETURN
421  ELSE IF( lquery ) THEN
422  RETURN
423  END IF
424 *
425 * Quick return if possible
426 *
427  m = 0
428  IF( n.EQ.0 ) THEN
429  RETURN
430  END IF
431 *
432  IF( n.EQ.1 ) THEN
433  IF( alleig .OR. indeig ) THEN
434  m = 1
435  w( 1 ) = real( a( 1, 1 ) )
436  ELSE IF( valeig ) THEN
437  IF( vl.LT.real( a( 1, 1 ) ) .AND. vu.GE.real( a( 1, 1 ) ) )
438  $ THEN
439  m = 1
440  w( 1 ) = real( a( 1, 1 ) )
441  END IF
442  END IF
443  IF( wantz )
444  $ z( 1, 1 ) = cone
445  RETURN
446  END IF
447 *
448 * Get machine constants.
449 *
450  safmin = slamch( 'Safe minimum' )
451  eps = slamch( 'Precision' )
452  smlnum = safmin / eps
453  bignum = one / smlnum
454  rmin = sqrt( smlnum )
455  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
456 *
457 * Scale matrix to allowable range, if necessary.
458 *
459  iscale = 0
460  abstll = abstol
461  IF( valeig ) THEN
462  vll = vl
463  vuu = vu
464  END IF
465  anrm = clanhe( 'M', uplo, n, a, lda, rwork )
466  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
467  iscale = 1
468  sigma = rmin / anrm
469  ELSE IF( anrm.GT.rmax ) THEN
470  iscale = 1
471  sigma = rmax / anrm
472  END IF
473  IF( iscale.EQ.1 ) THEN
474  IF( lower ) THEN
475  DO 10 j = 1, n
476  CALL csscal( n-j+1, sigma, a( j, j ), 1 )
477  10 CONTINUE
478  ELSE
479  DO 20 j = 1, n
480  CALL csscal( j, sigma, a( 1, j ), 1 )
481  20 CONTINUE
482  END IF
483  IF( abstol.GT.0 )
484  $ abstll = abstol*sigma
485  IF( valeig ) THEN
486  vll = vl*sigma
487  vuu = vu*sigma
488  END IF
489  END IF
490 *
491 * Call CHETRD_2STAGE to reduce Hermitian matrix to tridiagonal form.
492 *
493  indd = 1
494  inde = indd + n
495  indrwk = inde + n
496  indtau = 1
497  indhous = indtau + n
498  indwrk = indhous + lhtrd
499  llwork = lwork - indwrk + 1
500 *
501  CALL chetrd_2stage( jobz, uplo, n, a, lda, rwork( indd ),
502  $ rwork( inde ), work( indtau ),
503  $ work( indhous ), lhtrd, work( indwrk ),
504  $ llwork, iinfo )
505 *
506 * If all eigenvalues are desired and ABSTOL is less than or equal to
507 * zero, then call SSTERF or CUNGTR and CSTEQR. If this fails for
508 * some eigenvalue, then try SSTEBZ.
509 *
510  test = .false.
511  IF( indeig ) THEN
512  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
513  test = .true.
514  END IF
515  END IF
516  IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
517  CALL scopy( n, rwork( indd ), 1, w, 1 )
518  indee = indrwk + 2*n
519  IF( .NOT.wantz ) THEN
520  CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
521  CALL ssterf( n, w, rwork( indee ), info )
522  ELSE
523  CALL clacpy( 'A', n, n, a, lda, z, ldz )
524  CALL cungtr( uplo, n, z, ldz, work( indtau ),
525  $ work( indwrk ), llwork, iinfo )
526  CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
527  CALL csteqr( jobz, n, w, rwork( indee ), z, ldz,
528  $ rwork( indrwk ), info )
529  IF( info.EQ.0 ) THEN
530  DO 30 i = 1, n
531  ifail( i ) = 0
532  30 CONTINUE
533  END IF
534  END IF
535  IF( info.EQ.0 ) THEN
536  m = n
537  GO TO 40
538  END IF
539  info = 0
540  END IF
541 *
542 * Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
543 *
544  IF( wantz ) THEN
545  order = 'B'
546  ELSE
547  order = 'E'
548  END IF
549  indibl = 1
550  indisp = indibl + n
551  indiwk = indisp + n
552  CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
553  $ rwork( indd ), rwork( inde ), m, nsplit, w,
554  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
555  $ iwork( indiwk ), info )
556 *
557  IF( wantz ) THEN
558  CALL cstein( n, rwork( indd ), rwork( inde ), m, w,
559  $ iwork( indibl ), iwork( indisp ), z, ldz,
560  $ rwork( indrwk ), iwork( indiwk ), ifail, info )
561 *
562 * Apply unitary matrix used in reduction to tridiagonal
563 * form to eigenvectors returned by CSTEIN.
564 *
565  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
566  $ ldz, work( indwrk ), llwork, iinfo )
567  END IF
568 *
569 * If matrix was scaled, then rescale eigenvalues appropriately.
570 *
571  40 CONTINUE
572  IF( iscale.EQ.1 ) THEN
573  IF( info.EQ.0 ) THEN
574  imax = m
575  ELSE
576  imax = info - 1
577  END IF
578  CALL sscal( imax, one / sigma, w, 1 )
579  END IF
580 *
581 * If eigenvalues are not in order, then sort them, along with
582 * eigenvectors.
583 *
584  IF( wantz ) THEN
585  DO 60 j = 1, m - 1
586  i = 0
587  tmp1 = w( j )
588  DO 50 jj = j + 1, m
589  IF( w( jj ).LT.tmp1 ) THEN
590  i = jj
591  tmp1 = w( jj )
592  END IF
593  50 CONTINUE
594 *
595  IF( i.NE.0 ) THEN
596  itmp1 = iwork( indibl+i-1 )
597  w( i ) = w( j )
598  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
599  w( j ) = tmp1
600  iwork( indibl+j-1 ) = itmp1
601  CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
602  IF( info.NE.0 ) THEN
603  itmp1 = ifail( i )
604  ifail( i ) = ifail( j )
605  ifail( j ) = itmp1
606  END IF
607  END IF
608  60 CONTINUE
609  END IF
610 *
611 * Set WORK(1) to optimal complex workspace size.
612 *
613  work( 1 ) = lwmin
614 *
615  RETURN
616 *
617 * End of CHEEVX_2STAGE
618 *
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:149
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:86
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 csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:78
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
Definition: cswap.f:81
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,...
Definition: clanhe.f:124
subroutine chetrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
CHETRD_2STAGE
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:103
subroutine cunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMTR
Definition: cunmtr.f:172
subroutine cstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CSTEIN
Definition: cstein.f:182
subroutine csteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
CSTEQR
Definition: csteqr.f:132
subroutine cungtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
CUNGTR
Definition: cungtr.f:123
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:82
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:79
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:68
Here is the call graph for this function:
Here is the caller graph for this function: