LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cheevd ( character  JOBZ,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  W,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

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

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

Purpose:
 CHEEVD computes all eigenvalues and, optionally, eigenvectors of a
 complex Hermitian matrix A.  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.
[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 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is COMPLEX 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 REAL 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.
Date
November 2011
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 207 of file cheevd.f.

207 *
208 * -- LAPACK driver routine (version 3.4.0) --
209 * -- LAPACK is a software package provided by Univ. of Tennessee, --
210 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
211 * November 2011
212 *
213 * .. Scalar Arguments ..
214  CHARACTER jobz, uplo
215  INTEGER info, lda, liwork, lrwork, lwork, n
216 * ..
217 * .. Array Arguments ..
218  INTEGER iwork( * )
219  REAL rwork( * ), w( * )
220  COMPLEX a( lda, * ), work( * )
221 * ..
222 *
223 * =====================================================================
224 *
225 * .. Parameters ..
226  REAL zero, one
227  parameter ( zero = 0.0e0, one = 1.0e0 )
228  COMPLEX cone
229  parameter ( cone = ( 1.0e0, 0.0e0 ) )
230 * ..
231 * .. Local Scalars ..
232  LOGICAL lower, lquery, wantz
233  INTEGER iinfo, imax, inde, indrwk, indtau, indwk2,
234  $ indwrk, iscale, liopt, liwmin, llrwk, llwork,
235  $ llwrk2, lopt, lropt, lrwmin, lwmin
236  REAL anrm, bignum, eps, rmax, rmin, safmin, sigma,
237  $ smlnum
238 * ..
239 * .. External Functions ..
240  LOGICAL lsame
241  INTEGER ilaenv
242  REAL clanhe, slamch
243  EXTERNAL ilaenv, lsame, clanhe, slamch
244 * ..
245 * .. External Subroutines ..
246  EXTERNAL chetrd, clacpy, clascl, cstedc, cunmtr, sscal,
247  $ ssterf, xerbla
248 * ..
249 * .. Intrinsic Functions ..
250  INTRINSIC max, sqrt
251 * ..
252 * .. Executable Statements ..
253 *
254 * Test the input parameters.
255 *
256  wantz = lsame( jobz, 'V' )
257  lower = lsame( uplo, 'L' )
258  lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
259 *
260  info = 0
261  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
262  info = -1
263  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
264  info = -2
265  ELSE IF( n.LT.0 ) THEN
266  info = -3
267  ELSE IF( lda.LT.max( 1, n ) ) THEN
268  info = -5
269  END IF
270 *
271  IF( info.EQ.0 ) THEN
272  IF( n.LE.1 ) THEN
273  lwmin = 1
274  lrwmin = 1
275  liwmin = 1
276  lopt = lwmin
277  lropt = lrwmin
278  liopt = liwmin
279  ELSE
280  IF( wantz ) THEN
281  lwmin = 2*n + n*n
282  lrwmin = 1 + 5*n + 2*n**2
283  liwmin = 3 + 5*n
284  ELSE
285  lwmin = n + 1
286  lrwmin = n
287  liwmin = 1
288  END IF
289  lopt = max( lwmin, n +
290  $ ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 ) )
291  lropt = lrwmin
292  liopt = liwmin
293  END IF
294  work( 1 ) = lopt
295  rwork( 1 ) = lropt
296  iwork( 1 ) = liopt
297 *
298  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
299  info = -8
300  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
301  info = -10
302  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
303  info = -12
304  END IF
305  END IF
306 *
307  IF( info.NE.0 ) THEN
308  CALL xerbla( 'CHEEVD', -info )
309  RETURN
310  ELSE IF( lquery ) THEN
311  RETURN
312  END IF
313 *
314 * Quick return if possible
315 *
316  IF( n.EQ.0 )
317  $ RETURN
318 *
319  IF( n.EQ.1 ) THEN
320  w( 1 ) = a( 1, 1 )
321  IF( wantz )
322  $ a( 1, 1 ) = cone
323  RETURN
324  END IF
325 *
326 * Get machine constants.
327 *
328  safmin = slamch( 'Safe minimum' )
329  eps = slamch( 'Precision' )
330  smlnum = safmin / eps
331  bignum = one / smlnum
332  rmin = sqrt( smlnum )
333  rmax = sqrt( bignum )
334 *
335 * Scale matrix to allowable range, if necessary.
336 *
337  anrm = clanhe( 'M', uplo, n, a, lda, rwork )
338  iscale = 0
339  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
340  iscale = 1
341  sigma = rmin / anrm
342  ELSE IF( anrm.GT.rmax ) THEN
343  iscale = 1
344  sigma = rmax / anrm
345  END IF
346  IF( iscale.EQ.1 )
347  $ CALL clascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
348 *
349 * Call CHETRD to reduce Hermitian matrix to tridiagonal form.
350 *
351  inde = 1
352  indtau = 1
353  indwrk = indtau + n
354  indrwk = inde + n
355  indwk2 = indwrk + n*n
356  llwork = lwork - indwrk + 1
357  llwrk2 = lwork - indwk2 + 1
358  llrwk = lrwork - indrwk + 1
359  CALL chetrd( uplo, n, a, lda, w, rwork( inde ), work( indtau ),
360  $ work( indwrk ), llwork, iinfo )
361 *
362 * For eigenvalues only, call SSTERF. For eigenvectors, first call
363 * CSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
364 * tridiagonal matrix, then call CUNMTR to multiply it to the
365 * Householder transformations represented as Householder vectors in
366 * A.
367 *
368  IF( .NOT.wantz ) THEN
369  CALL ssterf( n, w, rwork( inde ), info )
370  ELSE
371  CALL cstedc( 'I', n, w, rwork( inde ), work( indwrk ), n,
372  $ work( indwk2 ), llwrk2, rwork( indrwk ), llrwk,
373  $ iwork, liwork, info )
374  CALL cunmtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
375  $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
376  CALL clacpy( 'A', n, n, work( indwrk ), n, a, lda )
377  END IF
378 *
379 * If matrix was scaled, then rescale eigenvalues appropriately.
380 *
381  IF( iscale.EQ.1 ) THEN
382  IF( info.EQ.0 ) THEN
383  imax = n
384  ELSE
385  imax = info - 1
386  END IF
387  CALL sscal( imax, one / sigma, w, 1 )
388  END IF
389 *
390  work( 1 ) = lopt
391  rwork( 1 ) = lropt
392  iwork( 1 ) = liopt
393 *
394  RETURN
395 *
396 * End of CHEEVD
397 *
subroutine clascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: clascl.f:145
subroutine cunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMTR
Definition: cunmtr.f:174
real function clanhe(NORM, UPLO, N, A, LDA, WORK)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix.
Definition: clanhe.f:126
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine chetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
CHETRD
Definition: chetrd.f:194
subroutine cstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CSTEDC
Definition: cstedc.f:214
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: