LAPACK 3.12.1
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.
!>
!> 
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 249 of file zhbevd_2stage.f.

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