LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zhbevx_2stage ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
integer  KD,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
complex*16, dimension( ldq, * )  Q,
integer  LDQ,
double precision  VL,
double precision  VU,
integer  IL,
integer  IU,
double precision  ABSTOL,
integer  M,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

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

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

Purpose:
 ZHBEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
 of a complex Hermitian band 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]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 0.
[in,out]AB
          AB is COMPLEX*16 array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix A, stored in the first KD+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).

          On exit, AB is overwritten by values generated during the
          reduction to tridiagonal form.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD + 1.
[out]Q
          Q is COMPLEX*16 array, dimension (LDQ, N)
          If JOBZ = 'V', the N-by-N unitary matrix used in the
                          reduction to tridiagonal form.
          If JOBZ = 'N', the array Q is not referenced.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  If JOBZ = 'V', then
          LDQ >= 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 AB 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)
          The first M elements contain the selected eigenvalues in
          ascending order.
[out]Z
          Z is COMPLEX*16 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*16 array, dimension (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, dimension) where
                                   dimension = (2KD+1)*N + KD*NTHREADS
                                   where KD is the size of the band.
                                   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 sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION 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
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 329 of file zhbevx_2stage.f.

329 *
330  IMPLICIT NONE
331 *
332 * -- LAPACK driver routine (version 3.7.0) --
333 * -- LAPACK is a software package provided by Univ. of Tennessee, --
334 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
335 * June 2016
336 *
337 * .. Scalar Arguments ..
338  CHARACTER jobz, range, uplo
339  INTEGER il, info, iu, kd, ldab, ldq, ldz, m, n, lwork
340  DOUBLE PRECISION abstol, vl, vu
341 * ..
342 * .. Array Arguments ..
343  INTEGER ifail( * ), iwork( * )
344  DOUBLE PRECISION rwork( * ), w( * )
345  COMPLEX*16 ab( ldab, * ), q( ldq, * ), work( * ),
346  $ z( ldz, * )
347 * ..
348 *
349 * =====================================================================
350 *
351 * .. Parameters ..
352  DOUBLE PRECISION zero, one
353  parameter ( zero = 0.0d0, one = 1.0d0 )
354  COMPLEX*16 czero, cone
355  parameter ( czero = ( 0.0d0, 0.0d0 ),
356  $ cone = ( 1.0d0, 0.0d0 ) )
357 * ..
358 * .. Local Scalars ..
359  LOGICAL alleig, indeig, lower, test, valeig, wantz,
360  $ lquery
361  CHARACTER order
362  INTEGER i, iinfo, imax, indd, inde, indee, indibl,
363  $ indisp, indiwk, indrwk, indwrk, iscale, itmp1,
364  $ llwork, lwmin, lhtrd, lwtrd, ib, indhous,
365  $ j, jj, nsplit
366  DOUBLE PRECISION abstll, anrm, bignum, eps, rmax, rmin, safmin,
367  $ sigma, smlnum, tmp1, vll, vuu
368  COMPLEX*16 ctmp1
369 * ..
370 * .. External Functions ..
371  LOGICAL lsame
372  INTEGER ilaenv
373  DOUBLE PRECISION dlamch, zlanhb
374  EXTERNAL lsame, dlamch, zlanhb, ilaenv
375 * ..
376 * .. External Subroutines ..
377  EXTERNAL dcopy, dscal, dstebz, dsterf, xerbla, zcopy,
379  $ zswap, zhetrd_hb2st
380 * ..
381 * .. Intrinsic Functions ..
382  INTRINSIC dble, max, min, sqrt
383 * ..
384 * .. Executable Statements ..
385 *
386 * Test the input parameters.
387 *
388  wantz = lsame( jobz, 'V' )
389  alleig = lsame( range, 'A' )
390  valeig = lsame( range, 'V' )
391  indeig = lsame( range, 'I' )
392  lower = lsame( uplo, 'L' )
393  lquery = ( lwork.EQ.-1 )
394 *
395  info = 0
396  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
397  info = -1
398  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
399  info = -2
400  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
401  info = -3
402  ELSE IF( n.LT.0 ) THEN
403  info = -4
404  ELSE IF( kd.LT.0 ) THEN
405  info = -5
406  ELSE IF( ldab.LT.kd+1 ) THEN
407  info = -7
408  ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
409  info = -9
410  ELSE
411  IF( valeig ) THEN
412  IF( n.GT.0 .AND. vu.LE.vl )
413  $ info = -11
414  ELSE IF( indeig ) THEN
415  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
416  info = -12
417  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
418  info = -13
419  END IF
420  END IF
421  END IF
422  IF( info.EQ.0 ) THEN
423  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
424  $ info = -18
425  END IF
426 *
427  IF( info.EQ.0 ) THEN
428  IF( n.LE.1 ) THEN
429  lwmin = 1
430  work( 1 ) = lwmin
431  ELSE
432  ib = ilaenv( 18, 'ZHETRD_HB2ST', jobz, n, kd, -1, -1 )
433  lhtrd = ilaenv( 19, 'ZHETRD_HB2ST', jobz, n, kd, ib, -1 )
434  lwtrd = ilaenv( 20, 'ZHETRD_HB2ST', jobz, n, kd, ib, -1 )
435  lwmin = lhtrd + lwtrd
436  work( 1 ) = lwmin
437  ENDIF
438 *
439  IF( lwork.LT.lwmin .AND. .NOT.lquery )
440  $ info = -20
441  END IF
442 *
443  IF( info.NE.0 ) THEN
444  CALL xerbla( 'ZHBEVX_2STAGE', -info )
445  RETURN
446  ELSE IF( lquery ) THEN
447  RETURN
448  END IF
449 *
450 * Quick return if possible
451 *
452  m = 0
453  IF( n.EQ.0 )
454  $ RETURN
455 *
456  IF( n.EQ.1 ) THEN
457  m = 1
458  IF( lower ) THEN
459  ctmp1 = ab( 1, 1 )
460  ELSE
461  ctmp1 = ab( kd+1, 1 )
462  END IF
463  tmp1 = dble( ctmp1 )
464  IF( valeig ) THEN
465  IF( .NOT.( vl.LT.tmp1 .AND. vu.GE.tmp1 ) )
466  $ m = 0
467  END IF
468  IF( m.EQ.1 ) THEN
469  w( 1 ) = dble( ctmp1 )
470  IF( wantz )
471  $ z( 1, 1 ) = cone
472  END IF
473  RETURN
474  END IF
475 *
476 * Get machine constants.
477 *
478  safmin = dlamch( 'Safe minimum' )
479  eps = dlamch( 'Precision' )
480  smlnum = safmin / eps
481  bignum = one / smlnum
482  rmin = sqrt( smlnum )
483  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
484 *
485 * Scale matrix to allowable range, if necessary.
486 *
487  iscale = 0
488  abstll = abstol
489  IF( valeig ) THEN
490  vll = vl
491  vuu = vu
492  ELSE
493  vll = zero
494  vuu = zero
495  END IF
496  anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
497  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
498  iscale = 1
499  sigma = rmin / anrm
500  ELSE IF( anrm.GT.rmax ) THEN
501  iscale = 1
502  sigma = rmax / anrm
503  END IF
504  IF( iscale.EQ.1 ) THEN
505  IF( lower ) THEN
506  CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
507  ELSE
508  CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
509  END IF
510  IF( abstol.GT.0 )
511  $ abstll = abstol*sigma
512  IF( valeig ) THEN
513  vll = vl*sigma
514  vuu = vu*sigma
515  END IF
516  END IF
517 *
518 * Call ZHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
519 *
520  indd = 1
521  inde = indd + n
522  indrwk = inde + n
523 *
524  indhous = 1
525  indwrk = indhous + lhtrd
526  llwork = lwork - indwrk + 1
527 *
528  CALL zhetrd_hb2st( 'N', jobz, uplo, n, kd, ab, ldab,
529  $ rwork( indd ), rwork( inde ), work( indhous ),
530  $ lhtrd, work( indwrk ), llwork, iinfo )
531 *
532 * If all eigenvalues are desired and ABSTOL is less than or equal
533 * to zero, then call DSTERF or ZSTEQR. If this fails for some
534 * eigenvalue, then try DSTEBZ.
535 *
536  test = .false.
537  IF (indeig) THEN
538  IF (il.EQ.1 .AND. iu.EQ.n) THEN
539  test = .true.
540  END IF
541  END IF
542  IF ((alleig .OR. test) .AND. (abstol.LE.zero)) THEN
543  CALL dcopy( n, rwork( indd ), 1, w, 1 )
544  indee = indrwk + 2*n
545  IF( .NOT.wantz ) THEN
546  CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
547  CALL dsterf( n, w, rwork( indee ), info )
548  ELSE
549  CALL zlacpy( 'A', n, n, q, ldq, z, ldz )
550  CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
551  CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
552  $ rwork( indrwk ), info )
553  IF( info.EQ.0 ) THEN
554  DO 10 i = 1, n
555  ifail( i ) = 0
556  10 CONTINUE
557  END IF
558  END IF
559  IF( info.EQ.0 ) THEN
560  m = n
561  GO TO 30
562  END IF
563  info = 0
564  END IF
565 *
566 * Otherwise, call DSTEBZ and, if eigenvectors are desired, ZSTEIN.
567 *
568  IF( wantz ) THEN
569  order = 'B'
570  ELSE
571  order = 'E'
572  END IF
573  indibl = 1
574  indisp = indibl + n
575  indiwk = indisp + n
576  CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
577  $ rwork( indd ), rwork( inde ), m, nsplit, w,
578  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
579  $ iwork( indiwk ), info )
580 *
581  IF( wantz ) THEN
582  CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
583  $ iwork( indibl ), iwork( indisp ), z, ldz,
584  $ rwork( indrwk ), iwork( indiwk ), ifail, info )
585 *
586 * Apply unitary matrix used in reduction to tridiagonal
587 * form to eigenvectors returned by ZSTEIN.
588 *
589  DO 20 j = 1, m
590  CALL zcopy( n, z( 1, j ), 1, work( 1 ), 1 )
591  CALL zgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
592  $ z( 1, j ), 1 )
593  20 CONTINUE
594  END IF
595 *
596 * If matrix was scaled, then rescale eigenvalues appropriately.
597 *
598  30 CONTINUE
599  IF( iscale.EQ.1 ) THEN
600  IF( info.EQ.0 ) THEN
601  imax = m
602  ELSE
603  imax = info - 1
604  END IF
605  CALL dscal( imax, one / sigma, w, 1 )
606  END IF
607 *
608 * If eigenvalues are not in order, then sort them, along with
609 * eigenvectors.
610 *
611  IF( wantz ) THEN
612  DO 50 j = 1, m - 1
613  i = 0
614  tmp1 = w( j )
615  DO 40 jj = j + 1, m
616  IF( w( jj ).LT.tmp1 ) THEN
617  i = jj
618  tmp1 = w( jj )
619  END IF
620  40 CONTINUE
621 *
622  IF( i.NE.0 ) THEN
623  itmp1 = iwork( indibl+i-1 )
624  w( i ) = w( j )
625  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
626  w( j ) = tmp1
627  iwork( indibl+j-1 ) = itmp1
628  CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
629  IF( info.NE.0 ) THEN
630  itmp1 = ifail( i )
631  ifail( i ) = ifail( j )
632  ifail( j ) = itmp1
633  END IF
634  END IF
635  50 CONTINUE
636  END IF
637 *
638 * Set WORK(1) to optimal workspace size.
639 *
640  work( 1 ) = lwmin
641 *
642  RETURN
643 *
644 * End of ZHBEVX_2STAGE
645 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: zlascl.f:145
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
Definition: zsteqr.f:134
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
Definition: zstein.f:184
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
double precision function zlanhb(NORM, UPLO, N, K, AB, LDAB, WORK)
ZLANHB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a Hermitian band matrix.
Definition: zlanhb.f:134

Here is the call graph for this function:

Here is the caller graph for this function: