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

◆ zhbevd_2stage()

subroutine zhbevd_2stage ( character  JOBZ,
character  UPLO,
integer  N,
integer  KD,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

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

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

Purpose:
 ZHBEVD_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
 a complex Hermitian band matrix A using the 2stage technique for
 the reduction to tridiagonal.  If eigenvectors are desired, it
 uses a divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
                  Not available in this release.
[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.  If UPLO = 'U', the first
          superdiagonal and the diagonal of the tridiagonal matrix T
          are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
          the diagonal and first subdiagonal of T are returned in the
          first two rows of AB.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD + 1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
          eigenvectors of the matrix A, with the i-th column of Z
          holding the eigenvector associated with W(i).
          If JOBZ = 'N', then Z is not referenced.
[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 (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, 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 (LRWORK)
          On exit, if INFO = 0, RWORK(1) returns the optimal LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of array RWORK.
          If N <= 1,               LRWORK must be at least 1.
          If JOBZ = 'N' and N > 1, LRWORK must be at least N.
          If JOBZ = 'V' and N > 1, LRWORK must be at least
                        1 + 5*N + 2*N**2.

          If LRWORK = -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]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ = 'N' or N <= 1, LIWORK must be at least 1.
          If JOBZ = 'V' and N > 1, LIWORK must be at least 3 + 5*N .

          If LIWORK = -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]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = i, the algorithm failed to converge; i
                off-diagonal elements of an intermediate tridiagonal
                form did not converge to zero.
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 257 of file zhbevd_2stage.f.

260*
261 IMPLICIT NONE
262*
263* -- LAPACK driver routine --
264* -- LAPACK is a software package provided by Univ. of Tennessee, --
265* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
266*
267* .. Scalar Arguments ..
268 CHARACTER JOBZ, UPLO
269 INTEGER INFO, KD, LDAB, LDZ, LIWORK, LRWORK, LWORK, N
270* ..
271* .. Array Arguments ..
272 INTEGER IWORK( * )
273 DOUBLE PRECISION RWORK( * ), W( * )
274 COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
275* ..
276*
277* =====================================================================
278*
279* .. Parameters ..
280 DOUBLE PRECISION ZERO, ONE
281 parameter( zero = 0.0d0, one = 1.0d0 )
282 COMPLEX*16 CZERO, CONE
283 parameter( czero = ( 0.0d0, 0.0d0 ),
284 $ cone = ( 1.0d0, 0.0d0 ) )
285* ..
286* .. Local Scalars ..
287 LOGICAL LOWER, LQUERY, WANTZ
288 INTEGER IINFO, IMAX, INDE, INDWK2, INDRWK, ISCALE,
289 $ LLWORK, INDWK, LHTRD, LWTRD, IB, INDHOUS,
290 $ LIWMIN, LLRWK, LLWK2, LRWMIN, LWMIN
291 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
292 $ SMLNUM
293* ..
294* .. External Functions ..
295 LOGICAL LSAME
296 INTEGER ILAENV2STAGE
297 DOUBLE PRECISION DLAMCH, ZLANHB
298 EXTERNAL lsame, dlamch, zlanhb, ilaenv2stage
299* ..
300* .. External Subroutines ..
301 EXTERNAL dscal, dsterf, xerbla, zgemm, zlacpy,
303* ..
304* .. Intrinsic Functions ..
305 INTRINSIC dble, sqrt
306* ..
307* .. Executable Statements ..
308*
309* Test the input parameters.
310*
311 wantz = lsame( jobz, 'V' )
312 lower = lsame( uplo, 'L' )
313 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 .OR. lrwork.EQ.-1 )
314*
315 info = 0
316 IF( n.LE.1 ) THEN
317 lwmin = 1
318 lrwmin = 1
319 liwmin = 1
320 ELSE
321 ib = ilaenv2stage( 2, 'ZHETRD_HB2ST', jobz, n, kd, -1, -1 )
322 lhtrd = ilaenv2stage( 3, 'ZHETRD_HB2ST', jobz, n, kd, ib, -1 )
323 lwtrd = ilaenv2stage( 4, 'ZHETRD_HB2ST', jobz, n, kd, ib, -1 )
324 IF( wantz ) THEN
325 lwmin = 2*n**2
326 lrwmin = 1 + 5*n + 2*n**2
327 liwmin = 3 + 5*n
328 ELSE
329 lwmin = max( n, lhtrd + lwtrd )
330 lrwmin = n
331 liwmin = 1
332 END IF
333 END IF
334 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
335 info = -1
336 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
337 info = -2
338 ELSE IF( n.LT.0 ) THEN
339 info = -3
340 ELSE IF( kd.LT.0 ) THEN
341 info = -4
342 ELSE IF( ldab.LT.kd+1 ) THEN
343 info = -6
344 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
345 info = -9
346 END IF
347*
348 IF( info.EQ.0 ) THEN
349 work( 1 ) = lwmin
350 rwork( 1 ) = lrwmin
351 iwork( 1 ) = liwmin
352*
353 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
354 info = -11
355 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
356 info = -13
357 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
358 info = -15
359 END IF
360 END IF
361*
362 IF( info.NE.0 ) THEN
363 CALL xerbla( 'ZHBEVD_2STAGE', -info )
364 RETURN
365 ELSE IF( lquery ) THEN
366 RETURN
367 END IF
368*
369* Quick return if possible
370*
371 IF( n.EQ.0 )
372 $ RETURN
373*
374 IF( n.EQ.1 ) THEN
375 w( 1 ) = dble( ab( 1, 1 ) )
376 IF( wantz )
377 $ z( 1, 1 ) = cone
378 RETURN
379 END IF
380*
381* Get machine constants.
382*
383 safmin = dlamch( 'Safe minimum' )
384 eps = dlamch( 'Precision' )
385 smlnum = safmin / eps
386 bignum = one / smlnum
387 rmin = sqrt( smlnum )
388 rmax = sqrt( bignum )
389*
390* Scale matrix to allowable range, if necessary.
391*
392 anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
393 iscale = 0
394 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
395 iscale = 1
396 sigma = rmin / anrm
397 ELSE IF( anrm.GT.rmax ) THEN
398 iscale = 1
399 sigma = rmax / anrm
400 END IF
401 IF( iscale.EQ.1 ) THEN
402 IF( lower ) THEN
403 CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
404 ELSE
405 CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
406 END IF
407 END IF
408*
409* Call ZHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
410*
411 inde = 1
412 indrwk = inde + n
413 llrwk = lrwork - indrwk + 1
414 indhous = 1
415 indwk = indhous + lhtrd
416 llwork = lwork - indwk + 1
417 indwk2 = indwk + n*n
418 llwk2 = lwork - indwk2 + 1
419*
420 CALL zhetrd_hb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
421 $ rwork( inde ), work( indhous ), lhtrd,
422 $ work( indwk ), llwork, iinfo )
423*
424* For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
425*
426 IF( .NOT.wantz ) THEN
427 CALL dsterf( n, w, rwork( inde ), info )
428 ELSE
429 CALL zstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
430 $ llwk2, rwork( indrwk ), llrwk, iwork, liwork,
431 $ info )
432 CALL zgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
433 $ work( indwk2 ), n )
434 CALL zlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
435 END IF
436*
437* If matrix was scaled, then rescale eigenvalues appropriately.
438*
439 IF( iscale.EQ.1 ) THEN
440 IF( info.EQ.0 ) THEN
441 imax = n
442 ELSE
443 imax = info - 1
444 END IF
445 CALL dscal( imax, one / sigma, w, 1 )
446 END IF
447*
448 work( 1 ) = lwmin
449 rwork( 1 ) = lrwmin
450 iwork( 1 ) = liwmin
451 RETURN
452*
453* End of ZHBEVD_2STAGE
454*
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
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 dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:86
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
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:143
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
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,...
Definition: zlanhb.f:132
subroutine zhetrd_hb2st(STAGE1, VECT, UPLO, N, KD, AB, LDAB, D, E, HOUS, LHOUS, WORK, LWORK, INFO)
ZHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
Definition: zhetrd_hb2st.F:230
subroutine zstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZSTEDC
Definition: zstedc.f:212
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:79
Here is the call graph for this function:
Here is the caller graph for this function: