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

◆ zhpev()

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.

Definition at line 136 of file zhpev.f.

138*
139* -- LAPACK driver routine --
140* -- LAPACK is a software package provided by Univ. of Tennessee, --
141* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
142*
143* .. Scalar Arguments ..
144 CHARACTER JOBZ, UPLO
145 INTEGER INFO, LDZ, N
146* ..
147* .. Array Arguments ..
148 DOUBLE PRECISION RWORK( * ), W( * )
149 COMPLEX*16 AP( * ), WORK( * ), Z( LDZ, * )
150* ..
151*
152* =====================================================================
153*
154* .. Parameters ..
155 DOUBLE PRECISION ZERO, ONE
156 parameter( zero = 0.0d0, one = 1.0d0 )
157* ..
158* .. Local Scalars ..
159 LOGICAL WANTZ
160 INTEGER IINFO, IMAX, INDE, INDRWK, INDTAU, INDWRK,
161 $ ISCALE
162 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
163 $ SMLNUM
164* ..
165* .. External Functions ..
166 LOGICAL LSAME
167 DOUBLE PRECISION DLAMCH, ZLANHP
168 EXTERNAL lsame, dlamch, zlanhp
169* ..
170* .. External Subroutines ..
171 EXTERNAL dscal, dsterf, xerbla, zdscal, zhptrd, zsteqr,
172 $ zupgtr
173* ..
174* .. Intrinsic Functions ..
175 INTRINSIC sqrt
176* ..
177* .. Executable Statements ..
178*
179* Test the input parameters.
180*
181 wantz = lsame( jobz, 'V' )
182*
183 info = 0
184 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
185 info = -1
186 ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
187 $ THEN
188 info = -2
189 ELSE IF( n.LT.0 ) THEN
190 info = -3
191 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
192 info = -7
193 END IF
194*
195 IF( info.NE.0 ) THEN
196 CALL xerbla( 'ZHPEV ', -info )
197 RETURN
198 END IF
199*
200* Quick return if possible
201*
202 IF( n.EQ.0 )
203 $ RETURN
204*
205 IF( n.EQ.1 ) THEN
206 w( 1 ) = dble( ap( 1 ) )
207 rwork( 1 ) = 1
208 IF( wantz )
209 $ z( 1, 1 ) = one
210 RETURN
211 END IF
212*
213* Get machine constants.
214*
215 safmin = dlamch( 'Safe minimum' )
216 eps = dlamch( 'Precision' )
217 smlnum = safmin / eps
218 bignum = one / smlnum
219 rmin = sqrt( smlnum )
220 rmax = sqrt( bignum )
221*
222* Scale matrix to allowable range, if necessary.
223*
224 anrm = zlanhp( 'M', uplo, n, ap, rwork )
225 iscale = 0
226 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
227 iscale = 1
228 sigma = rmin / anrm
229 ELSE IF( anrm.GT.rmax ) THEN
230 iscale = 1
231 sigma = rmax / anrm
232 END IF
233 IF( iscale.EQ.1 ) THEN
234 CALL zdscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
235 END IF
236*
237* Call ZHPTRD to reduce Hermitian packed matrix to tridiagonal form.
238*
239 inde = 1
240 indtau = 1
241 CALL zhptrd( uplo, n, ap, w, rwork( inde ), work( indtau ),
242 $ iinfo )
243*
244* For eigenvalues only, call DSTERF. For eigenvectors, first call
245* ZUPGTR to generate the orthogonal matrix, then call ZSTEQR.
246*
247 IF( .NOT.wantz ) THEN
248 CALL dsterf( n, w, rwork( inde ), info )
249 ELSE
250 indwrk = indtau + n
251 CALL zupgtr( uplo, n, ap, work( indtau ), z, ldz,
252 $ work( indwrk ), iinfo )
253 indrwk = inde + n
254 CALL zsteqr( jobz, n, w, rwork( inde ), z, ldz,
255 $ rwork( indrwk ), info )
256 END IF
257*
258* If matrix was scaled, then rescale eigenvalues appropriately.
259*
260 IF( iscale.EQ.1 ) THEN
261 IF( info.EQ.0 ) THEN
262 imax = n
263 ELSE
264 imax = info - 1
265 END IF
266 CALL dscal( imax, one / sigma, w, 1 )
267 END IF
268*
269 RETURN
270*
271* End of ZHPEV
272*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhptrd(uplo, n, ap, d, e, tau, info)
ZHPTRD
Definition zhptrd.f:151
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
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,...
Definition zlanhp.f:117
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:132
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine zupgtr(uplo, n, ap, tau, q, ldq, work, info)
ZUPGTR
Definition zupgtr.f:114
Here is the call graph for this function:
Here is the caller graph for this function: