LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhpev ( character  JOBZ,
character  UPLO,
integer  N,
complex*16, dimension( * )  AP,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  INFO 
)

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

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

Purpose:
 ZHPEV computes all the eigenvalues and, optionally, eigenvectors of a
 complex Hermitian matrix in packed storage.
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*16 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 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*16 array, dimension (max(1, 2*N-1))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(1, 3*N-2))
[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 140 of file zhpev.f.

140 *
141 * -- LAPACK driver routine (version 3.4.0) --
142 * -- LAPACK is a software package provided by Univ. of Tennessee, --
143 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144 * November 2011
145 *
146 * .. Scalar Arguments ..
147  CHARACTER jobz, uplo
148  INTEGER info, ldz, n
149 * ..
150 * .. Array Arguments ..
151  DOUBLE PRECISION rwork( * ), w( * )
152  COMPLEX*16 ap( * ), work( * ), z( ldz, * )
153 * ..
154 *
155 * =====================================================================
156 *
157 * .. Parameters ..
158  DOUBLE PRECISION zero, one
159  parameter ( zero = 0.0d0, one = 1.0d0 )
160 * ..
161 * .. Local Scalars ..
162  LOGICAL wantz
163  INTEGER iinfo, imax, inde, indrwk, indtau, indwrk,
164  $ iscale
165  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
166  $ smlnum
167 * ..
168 * .. External Functions ..
169  LOGICAL lsame
170  DOUBLE PRECISION dlamch, zlanhp
171  EXTERNAL lsame, dlamch, zlanhp
172 * ..
173 * .. External Subroutines ..
174  EXTERNAL dscal, dsterf, xerbla, zdscal, zhptrd, zsteqr,
175  $ zupgtr
176 * ..
177 * .. Intrinsic Functions ..
178  INTRINSIC sqrt
179 * ..
180 * .. Executable Statements ..
181 *
182 * Test the input parameters.
183 *
184  wantz = lsame( jobz, 'V' )
185 *
186  info = 0
187  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
188  info = -1
189  ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
190  $ THEN
191  info = -2
192  ELSE IF( n.LT.0 ) THEN
193  info = -3
194  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
195  info = -7
196  END IF
197 *
198  IF( info.NE.0 ) THEN
199  CALL xerbla( 'ZHPEV ', -info )
200  RETURN
201  END IF
202 *
203 * Quick return if possible
204 *
205  IF( n.EQ.0 )
206  $ RETURN
207 *
208  IF( n.EQ.1 ) THEN
209  w( 1 ) = ap( 1 )
210  rwork( 1 ) = 1
211  IF( wantz )
212  $ z( 1, 1 ) = one
213  RETURN
214  END IF
215 *
216 * Get machine constants.
217 *
218  safmin = dlamch( 'Safe minimum' )
219  eps = dlamch( 'Precision' )
220  smlnum = safmin / eps
221  bignum = one / smlnum
222  rmin = sqrt( smlnum )
223  rmax = sqrt( bignum )
224 *
225 * Scale matrix to allowable range, if necessary.
226 *
227  anrm = zlanhp( 'M', uplo, n, ap, rwork )
228  iscale = 0
229  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
230  iscale = 1
231  sigma = rmin / anrm
232  ELSE IF( anrm.GT.rmax ) THEN
233  iscale = 1
234  sigma = rmax / anrm
235  END IF
236  IF( iscale.EQ.1 ) THEN
237  CALL zdscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
238  END IF
239 *
240 * Call ZHPTRD to reduce Hermitian packed matrix to tridiagonal form.
241 *
242  inde = 1
243  indtau = 1
244  CALL zhptrd( uplo, n, ap, w, rwork( inde ), work( indtau ),
245  $ iinfo )
246 *
247 * For eigenvalues only, call DSTERF. For eigenvectors, first call
248 * ZUPGTR to generate the orthogonal matrix, then call ZSTEQR.
249 *
250  IF( .NOT.wantz ) THEN
251  CALL dsterf( n, w, rwork( inde ), info )
252  ELSE
253  indwrk = indtau + n
254  CALL zupgtr( uplo, n, ap, work( indtau ), z, ldz,
255  $ work( indwrk ), iinfo )
256  indrwk = inde + n
257  CALL zsteqr( jobz, n, w, rwork( inde ), z, ldz,
258  $ rwork( indrwk ), info )
259  END IF
260 *
261 * If matrix was scaled, then rescale eigenvalues appropriately.
262 *
263  IF( iscale.EQ.1 ) THEN
264  IF( info.EQ.0 ) THEN
265  imax = n
266  ELSE
267  imax = info - 1
268  END IF
269  CALL dscal( imax, one / sigma, w, 1 )
270  END IF
271 *
272  RETURN
273 *
274 * End of ZHPEV
275 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
double precision function zlanhp(NORM, UPLO, N, AP, WORK)
ZLANHP 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: zlanhp.f:119
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
Definition: zsteqr.f:134
subroutine zhptrd(UPLO, N, AP, D, E, TAU, INFO)
ZHPTRD
Definition: zhptrd.f:153
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zupgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
ZUPGTR
Definition: zupgtr.f:116

Here is the call graph for this function:

Here is the caller graph for this function: