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

◆ chpevd()

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.
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.

Definition at line 192 of file chpevd.f.

194*
195* -- LAPACK driver routine --
196* -- LAPACK is a software package provided by Univ. of Tennessee, --
197* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
198*
199* .. Scalar Arguments ..
200 CHARACTER JOBZ, UPLO
201 INTEGER INFO, LDZ, LIWORK, LRWORK, LWORK, N
202* ..
203* .. Array Arguments ..
204 INTEGER IWORK( * )
205 REAL RWORK( * ), W( * )
206 COMPLEX AP( * ), WORK( * ), Z( LDZ, * )
207* ..
208*
209* =====================================================================
210*
211* .. Parameters ..
212 REAL ZERO, ONE
213 parameter( zero = 0.0e+0, one = 1.0e+0 )
214 COMPLEX CONE
215 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
216* ..
217* .. Local Scalars ..
218 LOGICAL LQUERY, WANTZ
219 INTEGER IINFO, IMAX, INDE, INDRWK, INDTAU, INDWRK,
220 $ ISCALE, LIWMIN, LLRWK, LLWRK, LRWMIN, LWMIN
221 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
222 $ SMLNUM
223* ..
224* .. External Functions ..
225 LOGICAL LSAME
226 REAL CLANHP, SLAMCH, SROUNDUP_LWORK
228* ..
229* .. External Subroutines ..
230 EXTERNAL chptrd, csscal, cstedc, cupmtr, sscal, ssterf,
231 $ xerbla
232* ..
233* .. Intrinsic Functions ..
234 INTRINSIC sqrt
235* ..
236* .. Executable Statements ..
237*
238* Test the input parameters.
239*
240 wantz = lsame( jobz, 'V' )
241 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
242*
243 info = 0
244 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
245 info = -1
246 ELSE IF( .NOT.( lsame( uplo, 'L' ) .OR. lsame( uplo, 'U' ) ) )
247 $ THEN
248 info = -2
249 ELSE IF( n.LT.0 ) THEN
250 info = -3
251 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
252 info = -7
253 END IF
254*
255 IF( info.EQ.0 ) THEN
256 IF( n.LE.1 ) THEN
257 lwmin = 1
258 liwmin = 1
259 lrwmin = 1
260 ELSE
261 IF( wantz ) THEN
262 lwmin = 2*n
263 lrwmin = 1 + 5*n + 2*n**2
264 liwmin = 3 + 5*n
265 ELSE
266 lwmin = n
267 lrwmin = n
268 liwmin = 1
269 END IF
270 END IF
271 work( 1 ) = sroundup_lwork(lwmin)
272 rwork( 1 ) = lrwmin
273 iwork( 1 ) = liwmin
274*
275 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
276 info = -9
277 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
278 info = -11
279 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
280 info = -13
281 END IF
282 END IF
283*
284 IF( info.NE.0 ) THEN
285 CALL xerbla( 'CHPEVD', -info )
286 RETURN
287 ELSE IF( lquery ) THEN
288 RETURN
289 END IF
290*
291* Quick return if possible
292*
293 IF( n.EQ.0 )
294 $ RETURN
295*
296 IF( n.EQ.1 ) THEN
297 w( 1 ) = real( ap( 1 ) )
298 IF( wantz )
299 $ z( 1, 1 ) = cone
300 RETURN
301 END IF
302*
303* Get machine constants.
304*
305 safmin = slamch( 'Safe minimum' )
306 eps = slamch( 'Precision' )
307 smlnum = safmin / eps
308 bignum = one / smlnum
309 rmin = sqrt( smlnum )
310 rmax = sqrt( bignum )
311*
312* Scale matrix to allowable range, if necessary.
313*
314 anrm = clanhp( 'M', uplo, n, ap, rwork )
315 iscale = 0
316 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
317 iscale = 1
318 sigma = rmin / anrm
319 ELSE IF( anrm.GT.rmax ) THEN
320 iscale = 1
321 sigma = rmax / anrm
322 END IF
323 IF( iscale.EQ.1 ) THEN
324 CALL csscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
325 END IF
326*
327* Call CHPTRD to reduce Hermitian packed matrix to tridiagonal form.
328*
329 inde = 1
330 indtau = 1
331 indrwk = inde + n
332 indwrk = indtau + n
333 llwrk = lwork - indwrk + 1
334 llrwk = lrwork - indrwk + 1
335 CALL chptrd( uplo, n, ap, w, rwork( inde ), work( indtau ),
336 $ iinfo )
337*
338* For eigenvalues only, call SSTERF. For eigenvectors, first call
339* CUPGTR to generate the orthogonal matrix, then call CSTEDC.
340*
341 IF( .NOT.wantz ) THEN
342 CALL ssterf( n, w, rwork( inde ), info )
343 ELSE
344 CALL cstedc( 'I', n, w, rwork( inde ), z, ldz, work( indwrk ),
345 $ llwrk, rwork( indrwk ), llrwk, iwork, liwork,
346 $ info )
347 CALL cupmtr( 'L', uplo, 'N', n, n, ap, work( indtau ), z, ldz,
348 $ work( indwrk ), iinfo )
349 END IF
350*
351* If matrix was scaled, then rescale eigenvalues appropriately.
352*
353 IF( iscale.EQ.1 ) THEN
354 IF( info.EQ.0 ) THEN
355 imax = n
356 ELSE
357 imax = info - 1
358 END IF
359 CALL sscal( imax, one / sigma, w, 1 )
360 END IF
361*
362 work( 1 ) = sroundup_lwork(lwmin)
363 rwork( 1 ) = lrwmin
364 iwork( 1 ) = liwmin
365 RETURN
366*
367* End of CHPEVD
368*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chptrd(uplo, n, ap, d, e, tau, info)
CHPTRD
Definition chptrd.f:151
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clanhp(norm, uplo, n, ap, work)
CLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clanhp.f:117
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine cstedc(compz, n, d, e, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
CSTEDC
Definition cstedc.f:206
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
subroutine cupmtr(side, uplo, trans, m, n, ap, tau, c, ldc, work, info)
CUPMTR
Definition cupmtr.f:150
Here is the call graph for this function:
Here is the caller graph for this function: