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

◆ zheevd()

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

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

Purpose:
 ZHEEVD computes all eigenvalues and, optionally, eigenvectors of a
 complex Hermitian matrix A.  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.
[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 length of the array WORK.
          If N <= 1,                LWORK must be at least 1.
          If JOBZ  = 'N' and N > 1, LWORK must be at least N + 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

Definition at line 197 of file zheevd.f.

199*
200* -- LAPACK driver routine --
201* -- LAPACK is a software package provided by Univ. of Tennessee, --
202* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
203*
204* .. Scalar Arguments ..
205 CHARACTER JOBZ, UPLO
206 INTEGER INFO, LDA, LIWORK, LRWORK, LWORK, N
207* ..
208* .. Array Arguments ..
209 INTEGER IWORK( * )
210 DOUBLE PRECISION RWORK( * ), W( * )
211 COMPLEX*16 A( LDA, * ), WORK( * )
212* ..
213*
214* =====================================================================
215*
216* .. Parameters ..
217 DOUBLE PRECISION ZERO, ONE
218 parameter( zero = 0.0d0, one = 1.0d0 )
219 COMPLEX*16 CONE
220 parameter( cone = ( 1.0d0, 0.0d0 ) )
221* ..
222* .. Local Scalars ..
223 LOGICAL LOWER, LQUERY, WANTZ
224 INTEGER IINFO, IMAX, INDE, INDRWK, INDTAU, INDWK2,
225 $ INDWRK, ISCALE, LIOPT, LIWMIN, LLRWK, LLWORK,
226 $ LLWRK2, LOPT, LROPT, LRWMIN, LWMIN
227 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
228 $ SMLNUM
229* ..
230* .. External Functions ..
231 LOGICAL LSAME
232 INTEGER ILAENV
233 DOUBLE PRECISION DLAMCH, ZLANHE
234 EXTERNAL lsame, ilaenv, dlamch, zlanhe
235* ..
236* .. External Subroutines ..
237 EXTERNAL dscal, dsterf, xerbla, zhetrd, zlacpy, zlascl,
238 $ zstedc, zunmtr
239* ..
240* .. Intrinsic Functions ..
241 INTRINSIC max, sqrt
242* ..
243* .. Executable Statements ..
244*
245* Test the input parameters.
246*
247 wantz = lsame( jobz, 'V' )
248 lower = lsame( uplo, 'L' )
249 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
250*
251 info = 0
252 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
253 info = -1
254 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
255 info = -2
256 ELSE IF( n.LT.0 ) THEN
257 info = -3
258 ELSE IF( lda.LT.max( 1, n ) ) THEN
259 info = -5
260 END IF
261*
262 IF( info.EQ.0 ) THEN
263 IF( n.LE.1 ) THEN
264 lwmin = 1
265 lrwmin = 1
266 liwmin = 1
267 lopt = lwmin
268 lropt = lrwmin
269 liopt = liwmin
270 ELSE
271 IF( wantz ) THEN
272 lwmin = 2*n + n*n
273 lrwmin = 1 + 5*n + 2*n**2
274 liwmin = 3 + 5*n
275 ELSE
276 lwmin = n + 1
277 lrwmin = n
278 liwmin = 1
279 END IF
280 lopt = max( lwmin, n +
281 $ n*ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 ) )
282 lropt = lrwmin
283 liopt = liwmin
284 END IF
285 work( 1 ) = lopt
286 rwork( 1 ) = lropt
287 iwork( 1 ) = liopt
288*
289 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
290 info = -8
291 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
292 info = -10
293 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
294 info = -12
295 END IF
296 END IF
297*
298 IF( info.NE.0 ) THEN
299 CALL xerbla( 'ZHEEVD', -info )
300 RETURN
301 ELSE IF( lquery ) THEN
302 RETURN
303 END IF
304*
305* Quick return if possible
306*
307 IF( n.EQ.0 )
308 $ RETURN
309*
310 IF( n.EQ.1 ) THEN
311 w( 1 ) = dble( a( 1, 1 ) )
312 IF( wantz )
313 $ a( 1, 1 ) = cone
314 RETURN
315 END IF
316*
317* Get machine constants.
318*
319 safmin = dlamch( 'Safe minimum' )
320 eps = dlamch( 'Precision' )
321 smlnum = safmin / eps
322 bignum = one / smlnum
323 rmin = sqrt( smlnum )
324 rmax = sqrt( bignum )
325*
326* Scale matrix to allowable range, if necessary.
327*
328 anrm = zlanhe( 'M', uplo, n, a, lda, rwork )
329 iscale = 0
330 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
331 iscale = 1
332 sigma = rmin / anrm
333 ELSE IF( anrm.GT.rmax ) THEN
334 iscale = 1
335 sigma = rmax / anrm
336 END IF
337 IF( iscale.EQ.1 )
338 $ CALL zlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
339*
340* Call ZHETRD to reduce Hermitian matrix to tridiagonal form.
341*
342 inde = 1
343 indtau = 1
344 indwrk = indtau + n
345 indrwk = inde + n
346 indwk2 = indwrk + n*n
347 llwork = lwork - indwrk + 1
348 llwrk2 = lwork - indwk2 + 1
349 llrwk = lrwork - indrwk + 1
350 CALL zhetrd( uplo, n, a, lda, w, rwork( inde ), work( indtau ),
351 $ work( indwrk ), llwork, iinfo )
352*
353* For eigenvalues only, call DSTERF. For eigenvectors, first call
354* ZSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
355* tridiagonal matrix, then call ZUNMTR to multiply it to the
356* Householder transformations represented as Householder vectors in
357* A.
358*
359 IF( .NOT.wantz ) THEN
360 CALL dsterf( n, w, rwork( inde ), info )
361 ELSE
362 CALL zstedc( 'I', n, w, rwork( inde ), work( indwrk ), n,
363 $ work( indwk2 ), llwrk2, rwork( indrwk ), llrwk,
364 $ iwork, liwork, info )
365 CALL zunmtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
366 $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
367 CALL zlacpy( 'A', n, n, work( indwrk ), n, a, lda )
368 END IF
369*
370* If matrix was scaled, then rescale eigenvalues appropriately.
371*
372 IF( iscale.EQ.1 ) THEN
373 IF( info.EQ.0 ) THEN
374 imax = n
375 ELSE
376 imax = info - 1
377 END IF
378 CALL dscal( imax, one / sigma, w, 1 )
379 END IF
380*
381 work( 1 ) = lopt
382 rwork( 1 ) = lropt
383 iwork( 1 ) = liopt
384*
385 RETURN
386*
387* End of ZHEEVD
388*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhetrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
ZHETRD
Definition zhetrd.f:192
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
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 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:124
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
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:206
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine zunmtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
ZUNMTR
Definition zunmtr.f:171
Here is the call graph for this function:
Here is the caller graph for this function: