Visible to Intel only — GUID: GUID-D2916420-EF1C-4F16-B559-293A991EA303
Visible to Intel only — GUID: GUID-D2916420-EF1C-4F16-B559-293A991EA303
?hpevx
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix in packed storage.
Syntax
call chpevx(jobz, range, uplo, n, ap, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
call zhpevx(jobz, range, uplo, n, ap, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
call hpevx(ap, w [,uplo] [,z] [,vl] [,vu] [,il] [,iu] [,m] [,ifail] [,abstol] [,info])
Include Files
- mkl.fi, lapack.f90
Description
The routine computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A in packed storage. Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
Input Parameters
- jobz
-
CHARACTER*1. Must be 'N' or 'V'.
If job = 'N', then only eigenvalues are computed.
If job = 'V', then eigenvalues and eigenvectors are computed.
- range
-
CHARACTER*1. Must be 'A' or 'V' or 'I'.
If range = 'A', the routine computes all eigenvalues.
If range = 'V', the routine computes eigenvalues w(i) in the half-open interval: vl< w(i)≤vu.
If range = 'I', the routine computes eigenvalues with indices il to iu.
- uplo
-
CHARACTER*1. Must be 'U' or 'L'.
If uplo = 'U', ap stores the packed upper triangular part of A.
If uplo = 'L', ap stores the packed lower triangular part of A.
- n
-
INTEGER. The order of the matrix A (n≥ 0).
- ap, work
-
COMPLEX for chpevx
DOUBLE COMPLEX for zhpevx
Arrays:
Array ap(*) contains the packed upper or lower triangle of the Hermitian matrix A, as specified by uplo.
The size of ap must be at least max(1, n*(n+1)/2).
work(*) is a workspace array, size at least max(1, 2n).
- vl, vu
-
REAL for chpevx
DOUBLE PRECISION for zhpevx
If range = 'V', the lower and upper bounds of the interval to be searched for eigenvalues.
Constraint: vl< vu.
If range = 'A' or 'I', vl and vu are not referenced.
- il, iu
-
INTEGER.
If range = 'I', the indices in ascending order of the smallest and largest eigenvalues to be returned.
Constraint: 1 ≤il≤iu≤n, if n > 0; il=1 and iu=0 if n = 0.
If range = 'A' or 'V', il and iu are not referenced.
- abstol
-
REAL for chpevx
DOUBLE PRECISION for zhpevx
The absolute error tolerance to which each eigenvalue is required. See Application notes for details on error tolerance.
- ldz
-
INTEGER. The leading dimension of the output array z.
Constraints:
if jobz = 'N', then ldz≥ 1;
if jobz = 'V', then ldz≥ max(1, n).
- rwork
-
REAL for chpevx
DOUBLE PRECISION for zhpevx
Workspace array, size at least max(1, 7n).
- iwork
-
INTEGER. Workspace array, size at least max(1, 5n).
Output Parameters
- ap
-
On exit, this array is overwritten by the values generated during the reduction to tridiagonal form. The elements of the diagonal and the off-diagonal of the tridiagonal matrix overwrite the corresponding elements of A.
- m
-
INTEGER. The total number of eigenvalues found, 0 ≤m≤n.
0 ≤m≤n. If range = 'A', m = n, if range = 'I', m = iu-il+1, and if range = 'V' the exact value of m is not known in advance..
- w
-
REAL for chpevx
DOUBLE PRECISION for zhpevx
Array, size at least max(1, n).
If info = 0, contains the selected eigenvalues of the matrix A in ascending order.
- z
-
COMPLEX for chpevx
DOUBLE COMPLEX for zhpevx
Array z(ldz,*).
The second dimension of z must be at least max(1, m).
If jobz = 'V', then if info = 0, the first m columns of z contain the orthonormal eigenvectors of the matrix A corresponding to the selected eigenvalues, with the i-th column of z holding the eigenvector associated with w(i).
If an eigenvector fails to converge, then that column of z contains the latest approximation to the eigenvector, and the index of the eigenvector is returned in ifail.
If jobz = 'N', then z is not referenced.
Note: you must ensure that at least max(1,m) columns are supplied in the array z; if range = 'V', the exact value of m is not known in advance and an upper bound must be used.
- ifail
-
INTEGER.
Array, size at least max(1, n).
If jobz = 'V', then if info = 0, the first m elements of ifail are zero; if info > 0, the ifail contains the indices the eigenvectors that failed to converge.
If jobz = 'N', then ifail is not referenced.
- info
-
INTEGER.
If info = 0, the execution is successful.
If info = -i, the i-th parameter had an illegal value.
If info = i, then i eigenvectors failed to converge; their indices are stored in the array ifail.
LAPACK 95 Interface Notes
Routines in Fortran 95 interface have fewer arguments in the calling sequence than their FORTRAN 77 counterparts. For general conventions applied to skip redundant or restorable arguments, see LAPACK 95 Interface Conventions.
Specific details for the routine hpevx interface are the following:
- ap
-
Holds the array A of size (n*(n+1)/2).
- w
-
Holds the vector with the number of elements n.
- z
-
Holds the matrix Z of size (n, n), where the values n and m are significant.
- ifail
-
Holds the vector with the number of elements n.
- uplo
-
Must be 'U' or 'L'. The default value is 'U'.
- vl
-
Default value for this element is vl = -HUGE(vl).
- vu
-
Default value for this element is vu = HUGE(vl).
- il
-
Default value for this argument is il = 1.
- iu
-
Default value for this argument is iu = n.
- abstol
-
Default value for this element is abstol = 0.0_WP.
- jobz
-
Restored based on the presence of the argument z as follows:
jobz = 'V', if z is present,
jobz = 'N', if z is omitted
Note that there will be an error condition if ifail is present and z is omitted.
- range
-
Restored based on the presence of arguments vl, vu, il, iu as follows:
range = 'V', if one of or both vl and vu are present,
range = 'I', if one of or both il and iu are present,
range = 'A', if none of vl, vu, il, iu is present,
Note that there will be an error condition if one of or both vl and vu are present and at the same time one of or both il and iu are present.
Application Notes
An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a,b] of width less than or equal to abstol+ε*max(|a|,|b|), where ε is the machine precision.
If abstol is less than or equal to zero, then ε*||T||1 will be used in its place, where T is the tridiagonal matrix obtained by reducing A to tridiagonal form. Eigenvalues will be computed most accurately when abstol is set to twice the underflow threshold 2*?lamch('S'), not zero.
If this routine returns with info > 0, indicating that some eigenvectors did not converge, try setting abstol to 2*?lamch('S').