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

◆ zheevd_2stage()

subroutine zheevd_2stage ( character jobz,
character uplo,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
double precision, dimension( * ) w,
complex*16, dimension( * ) work,
integer lwork,
double precision, dimension( * ) rwork,
integer lrwork,
integer, dimension( * ) iwork,
integer liwork,
integer info )

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

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

Purpose:
!>
!> ZHEEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a
!> complex Hermitian 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,out]A
!>          A is COMPLEX*16 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, if JOBZ = 'V', then if INFO = 0, A contains the
!>          orthonormal eigenvectors of the matrix A.
!>          If JOBZ = 'N', then 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).
!> 
[out]W
!>          W is DOUBLE PRECISION array, dimension (N)
!>          If INFO = 0, the eigenvalues in ascending order.
!> 
[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 dimension of the array WORK.
!>          If N <= 1,               LWORK must be at least 1.
!>          If JOBZ = 'N' and N > 1, LWORK must be queried.
!>                                   LWORK = MAX(1, dimension) where
!>                                   dimension = max(stage1,stage2) + (KD+1)*N + N+1
!>                                             = N*KD + N*max(KD+1,FACTOPTNB)
!>                                               + max(2*KD*KD, KD*NTHREADS)
!>                                               + (KD+1)*N + N+1
!>                                   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 at least 2*N + N**2
!>
!>          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 the 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 the array IWORK.
!>          If N <= 1,                LIWORK must be at least 1.
!>          If JOBZ  = 'N' and 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 and JOBZ = 'N', then the algorithm failed
!>                to converge; i off-diagonal elements of an intermediate
!>                tridiagonal form did not converge to zero;
!>                if INFO = i and JOBZ = 'V', then the algorithm failed
!>                to compute an eigenvalue while working on the submatrix
!>                lying in rows and columns INFO/(N+1) through
!>                mod(INFO,N+1).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
Modified description of INFO. Sven, 16 Feb 05.
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA
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 243 of file zheevd_2stage.f.

246*
247 IMPLICIT NONE
248*
249* -- LAPACK driver routine --
250* -- LAPACK is a software package provided by Univ. of Tennessee, --
251* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
252*
253* .. Scalar Arguments ..
254 CHARACTER JOBZ, UPLO
255 INTEGER INFO, LDA, LIWORK, LRWORK, LWORK, N
256* ..
257* .. Array Arguments ..
258 INTEGER IWORK( * )
259 DOUBLE PRECISION RWORK( * ), W( * )
260 COMPLEX*16 A( LDA, * ), WORK( * )
261* ..
262*
263* =====================================================================
264*
265* .. Parameters ..
266 DOUBLE PRECISION ZERO, ONE
267 parameter( zero = 0.0d0, one = 1.0d0 )
268 COMPLEX*16 CONE
269 parameter( cone = ( 1.0d0, 0.0d0 ) )
270* ..
271* .. Local Scalars ..
272 LOGICAL LOWER, LQUERY, WANTZ
273 INTEGER IINFO, IMAX, INDE, INDRWK, INDTAU, INDWK2,
274 $ INDWRK, ISCALE, LIWMIN, LLRWK, LLWORK,
275 $ LLWRK2, LRWMIN, LWMIN,
276 $ LHTRD, LWTRD, KD, IB, INDHOUS
277
278
279 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
280 $ SMLNUM
281* ..
282* .. External Functions ..
283 LOGICAL LSAME
284 INTEGER ILAENV2STAGE
285 DOUBLE PRECISION DLAMCH, ZLANHE
286 EXTERNAL lsame, dlamch, zlanhe, ilaenv2stage
287* ..
288* .. External Subroutines ..
289 EXTERNAL dscal, dsterf, xerbla, zlacpy,
290 $ zlascl,
292* ..
293* .. Intrinsic Functions ..
294 INTRINSIC dble, max, sqrt
295* ..
296* .. Executable Statements ..
297*
298* Test the input parameters.
299*
300 wantz = lsame( jobz, 'V' )
301 lower = lsame( uplo, 'L' )
302 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
303*
304 info = 0
305 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
306 info = -1
307 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
308 info = -2
309 ELSE IF( n.LT.0 ) THEN
310 info = -3
311 ELSE IF( lda.LT.max( 1, n ) ) THEN
312 info = -5
313 END IF
314*
315 IF( info.EQ.0 ) THEN
316 IF( n.LE.1 ) THEN
317 lwmin = 1
318 lrwmin = 1
319 liwmin = 1
320 ELSE
321 kd = ilaenv2stage( 1, 'ZHETRD_2STAGE', jobz,
322 $ n, -1, -1, -1 )
323 ib = ilaenv2stage( 2, 'ZHETRD_2STAGE', jobz,
324 $ n, kd, -1, -1 )
325 lhtrd = ilaenv2stage( 3, 'ZHETRD_2STAGE', jobz,
326 $ n, kd, ib, -1 )
327 lwtrd = ilaenv2stage( 4, 'ZHETRD_2STAGE', jobz,
328 $ n, kd, ib, -1 )
329 IF( wantz ) THEN
330 lwmin = 2*n + n*n
331 lrwmin = 1 + 5*n + 2*n**2
332 liwmin = 3 + 5*n
333 ELSE
334 lwmin = n + 1 + lhtrd + lwtrd
335 lrwmin = n
336 liwmin = 1
337 END IF
338 END IF
339 work( 1 ) = lwmin
340 rwork( 1 ) = real( lrwmin )
341 iwork( 1 ) = liwmin
342*
343 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
344 info = -8
345 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
346 info = -10
347 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
348 info = -12
349 END IF
350 END IF
351*
352 IF( info.NE.0 ) THEN
353 CALL xerbla( 'ZHEEVD_2STAGE', -info )
354 RETURN
355 ELSE IF( lquery ) THEN
356 RETURN
357 END IF
358*
359* Quick return if possible
360*
361 IF( n.EQ.0 )
362 $ RETURN
363*
364 IF( n.EQ.1 ) THEN
365 w( 1 ) = dble( a( 1, 1 ) )
366 IF( wantz )
367 $ a( 1, 1 ) = cone
368 RETURN
369 END IF
370*
371* Get machine constants.
372*
373 safmin = dlamch( 'Safe minimum' )
374 eps = dlamch( 'Precision' )
375 smlnum = safmin / eps
376 bignum = one / smlnum
377 rmin = sqrt( smlnum )
378 rmax = sqrt( bignum )
379*
380* Scale matrix to allowable range, if necessary.
381*
382 anrm = zlanhe( 'M', uplo, n, a, lda, rwork )
383 iscale = 0
384 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
385 iscale = 1
386 sigma = rmin / anrm
387 ELSE IF( anrm.GT.rmax ) THEN
388 iscale = 1
389 sigma = rmax / anrm
390 END IF
391 IF( iscale.EQ.1 )
392 $ CALL zlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
393*
394* Call ZHETRD_2STAGE to reduce Hermitian matrix to tridiagonal form.
395*
396 inde = 1
397 indrwk = inde + n
398 llrwk = lrwork - indrwk + 1
399 indtau = 1
400 indhous = indtau + n
401 indwrk = indhous + lhtrd
402 llwork = lwork - indwrk + 1
403 indwk2 = indwrk + n*n
404 llwrk2 = lwork - indwk2 + 1
405*
406 CALL zhetrd_2stage( jobz, uplo, n, a, lda, w, rwork( inde ),
407 $ work( indtau ), work( indhous ), lhtrd,
408 $ work( indwrk ), llwork, iinfo )
409*
410* For eigenvalues only, call DSTERF. For eigenvectors, first call
411* ZSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
412* tridiagonal matrix, then call ZUNMTR to multiply it to the
413* Householder transformations represented as Householder vectors in
414* A.
415*
416 IF( .NOT.wantz ) THEN
417 CALL dsterf( n, w, rwork( inde ), info )
418 ELSE
419 CALL zstedc( 'I', n, w, rwork( inde ), work( indwrk ), n,
420 $ work( indwk2 ), llwrk2, rwork( indrwk ), llrwk,
421 $ iwork, liwork, info )
422 CALL zunmtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
423 $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
424 CALL zlacpy( 'A', n, n, work( indwrk ), n, a, lda )
425 END IF
426*
427* If matrix was scaled, then rescale eigenvalues appropriately.
428*
429 IF( iscale.EQ.1 ) THEN
430 IF( info.EQ.0 ) THEN
431 imax = n
432 ELSE
433 imax = info - 1
434 END IF
435 CALL dscal( imax, one / sigma, w, 1 )
436 END IF
437*
438 work( 1 ) = lwmin
439 rwork( 1 ) = real( lrwmin )
440 iwork( 1 ) = liwmin
441*
442 RETURN
443*
444* End of ZHEEVD_2STAGE
445*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhetrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
ZHETRD_2STAGE
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 zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:122
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
subroutine zunmtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
ZUNMTR
Definition zunmtr.f:170
Here is the call graph for this function:
Here is the caller graph for this function: