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

◆ zhpevd()

subroutine zhpevd ( 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,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

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

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

Purpose:
 ZHPEVD 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*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,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 DOUBLE PRECISION 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.

Definition at line 198 of file zhpevd.f.

200*
201* -- LAPACK driver routine --
202* -- LAPACK is a software package provided by Univ. of Tennessee, --
203* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
204*
205* .. Scalar Arguments ..
206 CHARACTER JOBZ, UPLO
207 INTEGER INFO, LDZ, LIWORK, LRWORK, LWORK, N
208* ..
209* .. Array Arguments ..
210 INTEGER IWORK( * )
211 DOUBLE PRECISION RWORK( * ), W( * )
212 COMPLEX*16 AP( * ), WORK( * ), Z( LDZ, * )
213* ..
214*
215* =====================================================================
216*
217* .. Parameters ..
218 DOUBLE PRECISION ZERO, ONE
219 parameter( zero = 0.0d+0, one = 1.0d+0 )
220 COMPLEX*16 CONE
221 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
222* ..
223* .. Local Scalars ..
224 LOGICAL LQUERY, WANTZ
225 INTEGER IINFO, IMAX, INDE, INDRWK, INDTAU, INDWRK,
226 $ ISCALE, LIWMIN, LLRWK, LLWRK, LRWMIN, LWMIN
227 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
228 $ SMLNUM
229* ..
230* .. External Functions ..
231 LOGICAL LSAME
232 DOUBLE PRECISION DLAMCH, ZLANHP
233 EXTERNAL lsame, dlamch, zlanhp
234* ..
235* .. External Subroutines ..
236 EXTERNAL dscal, dsterf, xerbla, zdscal, zhptrd, zstedc,
237 $ zupmtr
238* ..
239* .. Intrinsic Functions ..
240 INTRINSIC sqrt
241* ..
242* .. Executable Statements ..
243*
244* Test the input parameters.
245*
246 wantz = lsame( jobz, 'V' )
247 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
248*
249 info = 0
250 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
251 info = -1
252 ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
253 $ THEN
254 info = -2
255 ELSE IF( n.LT.0 ) THEN
256 info = -3
257 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
258 info = -7
259 END IF
260*
261 IF( info.EQ.0 ) THEN
262 IF( n.LE.1 ) THEN
263 lwmin = 1
264 liwmin = 1
265 lrwmin = 1
266 ELSE
267 IF( wantz ) THEN
268 lwmin = 2*n
269 lrwmin = 1 + 5*n + 2*n**2
270 liwmin = 3 + 5*n
271 ELSE
272 lwmin = n
273 lrwmin = n
274 liwmin = 1
275 END IF
276 END IF
277 work( 1 ) = lwmin
278 rwork( 1 ) = lrwmin
279 iwork( 1 ) = liwmin
280*
281 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
282 info = -9
283 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
284 info = -11
285 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
286 info = -13
287 END IF
288 END IF
289*
290 IF( info.NE.0 ) THEN
291 CALL xerbla( 'ZHPEVD', -info )
292 RETURN
293 ELSE IF( lquery ) THEN
294 RETURN
295 END IF
296*
297* Quick return if possible
298*
299 IF( n.EQ.0 )
300 $ RETURN
301*
302 IF( n.EQ.1 ) THEN
303 w( 1 ) = dble( ap( 1 ) )
304 IF( wantz )
305 $ z( 1, 1 ) = cone
306 RETURN
307 END IF
308*
309* Get machine constants.
310*
311 safmin = dlamch( 'Safe minimum' )
312 eps = dlamch( 'Precision' )
313 smlnum = safmin / eps
314 bignum = one / smlnum
315 rmin = sqrt( smlnum )
316 rmax = sqrt( bignum )
317*
318* Scale matrix to allowable range, if necessary.
319*
320 anrm = zlanhp( 'M', uplo, n, ap, rwork )
321 iscale = 0
322 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
323 iscale = 1
324 sigma = rmin / anrm
325 ELSE IF( anrm.GT.rmax ) THEN
326 iscale = 1
327 sigma = rmax / anrm
328 END IF
329 IF( iscale.EQ.1 ) THEN
330 CALL zdscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
331 END IF
332*
333* Call ZHPTRD to reduce Hermitian packed matrix to tridiagonal form.
334*
335 inde = 1
336 indtau = 1
337 indrwk = inde + n
338 indwrk = indtau + n
339 llwrk = lwork - indwrk + 1
340 llrwk = lrwork - indrwk + 1
341 CALL zhptrd( uplo, n, ap, w, rwork( inde ), work( indtau ),
342 $ iinfo )
343*
344* For eigenvalues only, call DSTERF. For eigenvectors, first call
345* ZUPGTR to generate the orthogonal matrix, then call ZSTEDC.
346*
347 IF( .NOT.wantz ) THEN
348 CALL dsterf( n, w, rwork( inde ), info )
349 ELSE
350 CALL zstedc( 'I', n, w, rwork( inde ), z, ldz, work( indwrk ),
351 $ llwrk, rwork( indrwk ), llrwk, iwork, liwork,
352 $ info )
353 CALL zupmtr( 'L', uplo, 'N', n, n, ap, work( indtau ), z, ldz,
354 $ work( indwrk ), iinfo )
355 END IF
356*
357* If matrix was scaled, then rescale eigenvalues appropriately.
358*
359 IF( iscale.EQ.1 ) THEN
360 IF( info.EQ.0 ) THEN
361 imax = n
362 ELSE
363 imax = info - 1
364 END IF
365 CALL dscal( imax, one / sigma, w, 1 )
366 END IF
367*
368 work( 1 ) = lwmin
369 rwork( 1 ) = lrwmin
370 iwork( 1 ) = liwmin
371 RETURN
372*
373* End of ZHPEVD
374*
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:86
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:78
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
subroutine zupmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
ZUPMTR
Definition: zupmtr.f:150
subroutine zhptrd(UPLO, N, AP, D, E, TAU, INFO)
ZHPTRD
Definition: zhptrd.f:151
subroutine zstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZSTEDC
Definition: zstedc.f:212
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:79
Here is the call graph for this function:
Here is the caller graph for this function: