LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine chpevd ( character  JOBZ,
character  UPLO,
integer  N,
complex, dimension( * )  AP,
real, dimension( * )  W,
complex, dimension( ldz, * )  Z,
integer  LDZ,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

CHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

Purpose:
 CHPEVD computes all the eigenvalues and, optionally, eigenvectors of
 a complex Hermitian matrix A in packed storage.  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]AP
          AP is COMPLEX array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, AP is overwritten by values generated during the
          reduction to tridiagonal form.  If UPLO = 'U', the diagonal
          and first superdiagonal of the tridiagonal matrix T overwrite
          the corresponding elements of A, and if UPLO = 'L', the
          diagonal and first subdiagonal of T overwrite the
          corresponding elements of A.
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX 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 array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of array WORK.
          If N <= 1,               LWORK must be at least 1.
          If JOBZ = 'N' and N > 1, LWORK must be at least N.
          If JOBZ = 'V' and N > 1, LWORK must be at least 2*N.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required 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 (MAX(1,LRWORK))
          On exit, if INFO = 0, RWORK(1) returns the required 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 required 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 required 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 required 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.
Date
November 2011

Definition at line 202 of file chpevd.f.

202 *
203 * -- LAPACK driver routine (version 3.4.0) --
204 * -- LAPACK is a software package provided by Univ. of Tennessee, --
205 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
206 * November 2011
207 *
208 * .. Scalar Arguments ..
209  CHARACTER jobz, uplo
210  INTEGER info, ldz, liwork, lrwork, lwork, n
211 * ..
212 * .. Array Arguments ..
213  INTEGER iwork( * )
214  REAL rwork( * ), w( * )
215  COMPLEX ap( * ), work( * ), z( ldz, * )
216 * ..
217 *
218 * =====================================================================
219 *
220 * .. Parameters ..
221  REAL zero, one
222  parameter ( zero = 0.0e+0, one = 1.0e+0 )
223  COMPLEX cone
224  parameter ( cone = ( 1.0e+0, 0.0e+0 ) )
225 * ..
226 * .. Local Scalars ..
227  LOGICAL lquery, wantz
228  INTEGER iinfo, imax, inde, indrwk, indtau, indwrk,
229  $ iscale, liwmin, llrwk, llwrk, lrwmin, lwmin
230  REAL anrm, bignum, eps, rmax, rmin, safmin, sigma,
231  $ smlnum
232 * ..
233 * .. External Functions ..
234  LOGICAL lsame
235  REAL clanhp, slamch
236  EXTERNAL lsame, clanhp, slamch
237 * ..
238 * .. External Subroutines ..
239  EXTERNAL chptrd, csscal, cstedc, cupmtr, sscal, ssterf,
240  $ xerbla
241 * ..
242 * .. Intrinsic Functions ..
243  INTRINSIC sqrt
244 * ..
245 * .. Executable Statements ..
246 *
247 * Test the input parameters.
248 *
249  wantz = lsame( jobz, 'V' )
250  lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
251 *
252  info = 0
253  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
254  info = -1
255  ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
256  $ THEN
257  info = -2
258  ELSE IF( n.LT.0 ) THEN
259  info = -3
260  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
261  info = -7
262  END IF
263 *
264  IF( info.EQ.0 ) THEN
265  IF( n.LE.1 ) THEN
266  lwmin = 1
267  liwmin = 1
268  lrwmin = 1
269  ELSE
270  IF( wantz ) THEN
271  lwmin = 2*n
272  lrwmin = 1 + 5*n + 2*n**2
273  liwmin = 3 + 5*n
274  ELSE
275  lwmin = n
276  lrwmin = n
277  liwmin = 1
278  END IF
279  END IF
280  work( 1 ) = lwmin
281  rwork( 1 ) = lrwmin
282  iwork( 1 ) = liwmin
283 *
284  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
285  info = -9
286  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
287  info = -11
288  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
289  info = -13
290  END IF
291  END IF
292 *
293  IF( info.NE.0 ) THEN
294  CALL xerbla( 'CHPEVD', -info )
295  RETURN
296  ELSE IF( lquery ) THEN
297  RETURN
298  END IF
299 *
300 * Quick return if possible
301 *
302  IF( n.EQ.0 )
303  $ RETURN
304 *
305  IF( n.EQ.1 ) THEN
306  w( 1 ) = ap( 1 )
307  IF( wantz )
308  $ z( 1, 1 ) = cone
309  RETURN
310  END IF
311 *
312 * Get machine constants.
313 *
314  safmin = slamch( 'Safe minimum' )
315  eps = slamch( 'Precision' )
316  smlnum = safmin / eps
317  bignum = one / smlnum
318  rmin = sqrt( smlnum )
319  rmax = sqrt( bignum )
320 *
321 * Scale matrix to allowable range, if necessary.
322 *
323  anrm = clanhp( 'M', uplo, n, ap, rwork )
324  iscale = 0
325  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
326  iscale = 1
327  sigma = rmin / anrm
328  ELSE IF( anrm.GT.rmax ) THEN
329  iscale = 1
330  sigma = rmax / anrm
331  END IF
332  IF( iscale.EQ.1 ) THEN
333  CALL csscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
334  END IF
335 *
336 * Call CHPTRD to reduce Hermitian packed matrix to tridiagonal form.
337 *
338  inde = 1
339  indtau = 1
340  indrwk = inde + n
341  indwrk = indtau + n
342  llwrk = lwork - indwrk + 1
343  llrwk = lrwork - indrwk + 1
344  CALL chptrd( uplo, n, ap, w, rwork( inde ), work( indtau ),
345  $ iinfo )
346 *
347 * For eigenvalues only, call SSTERF. For eigenvectors, first call
348 * CUPGTR to generate the orthogonal matrix, then call CSTEDC.
349 *
350  IF( .NOT.wantz ) THEN
351  CALL ssterf( n, w, rwork( inde ), info )
352  ELSE
353  CALL cstedc( 'I', n, w, rwork( inde ), z, ldz, work( indwrk ),
354  $ llwrk, rwork( indrwk ), llrwk, iwork, liwork,
355  $ info )
356  CALL cupmtr( 'L', uplo, 'N', n, n, ap, work( indtau ), z, ldz,
357  $ work( indwrk ), iinfo )
358  END IF
359 *
360 * If matrix was scaled, then rescale eigenvalues appropriately.
361 *
362  IF( iscale.EQ.1 ) THEN
363  IF( info.EQ.0 ) THEN
364  imax = n
365  ELSE
366  imax = info - 1
367  END IF
368  CALL sscal( imax, one / sigma, w, 1 )
369  END IF
370 *
371  work( 1 ) = lwmin
372  rwork( 1 ) = lrwmin
373  iwork( 1 ) = liwmin
374  RETURN
375 *
376 * End of CHPEVD
377 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CSTEDC
Definition: cstedc.f:214
subroutine chptrd(UPLO, N, AP, D, E, TAU, INFO)
CHPTRD
Definition: chptrd.f:153
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
real function clanhp(NORM, UPLO, N, AP, WORK)
CLANHP 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 supplied in packed form.
Definition: clanhp.f:119
subroutine cupmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
CUPMTR
Definition: cupmtr.f:152
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
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: