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

◆ cheevd()

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.
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 203 of file cheevd.f.

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