LAPACK++
2022.05.00
LAPACK C++ API
|
Functions | |
int64_t | lapack::heev (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double *W) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A. More... | |
int64_t | lapack::heev (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float *W) |
int64_t | lapack::heev_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double *W) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal. More... | |
int64_t | lapack::heev_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float *W) |
int64_t | lapack::heevd (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double *W) |
Computes all eigenvalues and, optionally, eigenvectors of a. More... | |
int64_t | lapack::heevd (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float *W) |
int64_t | lapack::heevd_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double *W) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal. More... | |
int64_t | lapack::heevd_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float *W) |
int64_t | lapack::heevr (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *nfound, double *W, std::complex< double > *Z, int64_t ldz, int64_t *isuppz) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A. More... | |
int64_t | lapack::heevr (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, std::complex< float > *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::heevr_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *nfound, double *W, std::complex< double > *Z, int64_t ldz, int64_t *isuppz) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal. More... | |
int64_t | lapack::heevr_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, std::complex< float > *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::heevx (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *nfound, double *W, std::complex< double > *Z, int64_t ldz, int64_t *ifail) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A. More... | |
int64_t | lapack::heevx (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, std::complex< float > *Z, int64_t ldz, int64_t *ifail) |
int64_t | lapack::heevx_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *nfound, double *W, std::complex< double > *Z, int64_t ldz, int64_t *ifail) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal. More... | |
int64_t | lapack::heevx_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, std::complex< float > *Z, int64_t ldz, int64_t *ifail) |
int64_t | lapack::syev (lapack::Job jobz, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double *W) |
int64_t | lapack::syev (lapack::Job jobz, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float *W) |
int64_t | lapack::syev_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double *W) |
int64_t | lapack::syev_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float *W) |
int64_t | lapack::syevd (lapack::Job jobz, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double *W) |
int64_t | lapack::syevd (lapack::Job jobz, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float *W) |
int64_t | lapack::syevd_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double *W) |
int64_t | lapack::syevd_2stage (lapack::Job jobz, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float *W) |
int64_t | lapack::syevr (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *m, double *W, double *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::syevr (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, float *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::syevr_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *m, double *W, double *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::syevr_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, float *Z, int64_t ldz, int64_t *isuppz) |
int64_t | lapack::syevx (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *m, double *W, double *Z, int64_t ldz, int64_t *ifail) |
int64_t | lapack::syevx (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, float *Z, int64_t ldz, int64_t *ifail) |
int64_t | lapack::syevx_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double vl, double vu, int64_t il, int64_t iu, double abstol, int64_t *m, double *W, double *Z, int64_t ldz, int64_t *ifail) |
int64_t | lapack::syevx_2stage (lapack::Job jobz, lapack::Range range, lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float vl, float vu, int64_t il, int64_t iu, float abstol, int64_t *m, float *W, float *Z, int64_t ldz, int64_t *ifail) |
int64_t lapack::heev | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syev
.
[in] | jobz |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | W | The vector W of length n. If successful, the eigenvalues in ascending order. |
int64_t lapack::heev_2stage | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syev_2stage
.
[in] | jobz |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | W | The vector W of length n. If successful, the eigenvalues in ascending order. |
All details about the 2-stage techniques are available in:
Azzam Haidar, Hatem Ltaief, and Jack Dongarra. Parallel reduction to condensed forms for symmetric eigenvalue problems using aggregated fine-grained and memory-aware kernels. In Proceedings of 2011 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '11), New York, NY, USA, Article 8, 11 pages. http://doi.acm.org/10.1145/2063384.2063394
A. Haidar, J. Kurzak, P. Luszczek, 2013. An improved parallel singular value algorithm and its implementation for multicore hardware, In Proceedings of 2013 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '13). Denver, Colorado, USA, 2013. Article 90, 12 pages. http://doi.acm.org/10.1145/2503210.2503292
A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. A novel hybrid CPU-GPU generalized eigensolver for electronic structure calculations based on fine-grained memory aware tasks. International Journal of High Performance Computing Applications. Volume 28 Issue 2, Pages 196-209, May 2014. http://hpc.sagepub.com/content/28/2/196
int64_t lapack::heevd | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
Computes all eigenvalues and, optionally, eigenvectors of 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.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syevd
.
[in] | jobz |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | W | The vector W of length n. If successful, the eigenvalues in ascending order. |
int64_t lapack::heevd_2stage | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal.
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.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syevd_2stage
.
[in] | jobz |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | W | The vector W of length n. If successful, the eigenvalues in ascending order. |
All details about the 2-stage techniques are available in:
Azzam Haidar, Hatem Ltaief, and Jack Dongarra. Parallel reduction to condensed forms for symmetric eigenvalue problems using aggregated fine-grained and memory-aware kernels. In Proceedings of 2011 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '11), New York, NY, USA, Article 8, 11 pages. http://doi.acm.org/10.1145/2063384.2063394
A. Haidar, J. Kurzak, P. Luszczek, 2013. An improved parallel singular value algorithm and its implementation for multicore hardware, In Proceedings of 2013 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '13). Denver, Colorado, USA, 2013. Article 90, 12 pages. http://doi.acm.org/10.1145/2503210.2503292
A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. A novel hybrid CPU-GPU generalized eigensolver for electronic structure calculations based on fine-grained memory aware tasks. International Journal of High Performance Computing Applications. Volume 28 Issue 2, Pages 196-209, May 2014. http://hpc.sagepub.com/content/28/2/196
int64_t lapack::heevr | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
std::complex< double > * | Z, | ||
int64_t | ldz, | ||
int64_t * | isuppz | ||
) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A.
Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
heevr
first reduces the matrix A to tridiagonal form T with a call to lapack::hetrd
. Then, whenever possible, heevr
calls lapack::stemr
to compute eigenspectrum using Relatively Robust Representations. lapack::stemr
computes eigenvalues by the dqds algorithm, while orthogonal eigenvectors are computed from various "good" \(L D L^T\) representations (also known as Relatively Robust Representations). Gram-Schmidt orthogonalization is avoided as far as possible. More specifically, the various steps of the algorithm are as follows.
For each unreduced block (submatrix) of T,
(a) Compute \(T - \sigma I = L D L^T\), so that L and D define all the wanted eigenvalues to high relative accuracy. This means that small relative changes in the entries of D and L cause only small relative changes in the eigenvalues and eigenvectors. The standard (unfactored) representation of the tridiagonal matrix T does not have this property in general.
(b) Compute the eigenvalues to suitable accuracy. If the eigenvectors are desired, the algorithm attains full accuracy of the computed eigenvalues only right before the corresponding vectors have to be computed, see steps c) and d).
(c) For each cluster of close eigenvalues, select a new shift close to the cluster, find a new factorization, and refine the shifted eigenvalues to suitable accuracy.
(d) For each eigenvalue with a large enough relative separation compute the corresponding eigenvector by forming a rank revealing twisted factorization. Go back to (c) for any clusters that remain.
The desired accuracy of the output can be specified by the input parameter abstol.
For more details, see lapack::stemr
documentation and:
Note 1 : heevr
calls lapack::stemr
when the full spectrum is requested on machines which conform to the ieee-754 floating point standard. heevr
calls lapack::stebz
and lapack::stein
on non-IEEE machines and when partial spectrum requests are made.
Normal execution of lapack::stemr
may create NaNs and infinities and hence may abort due to a floating point exception in environments which do not handle NaNs and infinities in the IEEE standard default manner.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syevr
.
[in] | jobz |
|
[in] | range |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | vl | If range=Value, the lower bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | vu | If range=Value, the upper bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | il | If range=Index, the index of the smallest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | iu | If range=Index, the index of the largest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | abstol | The absolute error tolerance for the eigenvalues. 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 + eps * max( |a|,|b| ), where eps is the machine precision. If abstol is less than or equal to zero, then eps*|T| will be used in its place, where |T| is the 1-norm of the tridiagonal matrix obtained by reducing A to tridiagonal form. See "Computing Small Singular Values of Bidiagonal Matrices with Guaranteed High Relative Accuracy," by Demmel and Kahan, LAPACK Working Note #3. If high relative accuracy is important, set abstol to DLAMCH( 'Safe minimum' ). Doing so will guarantee that eigenvalues are computed to high relative accuracy when possible in future releases. The current code does not make any guarantees about high relative accuracy, but future releases will. See J. Barlow and J. Demmel, "Computing Accurate Eigensystems of Scaled Diagonally Dominant Matrices", LAPACK Working Note #7, for a discussion of which matrices define their eigenvalues to high relative accuracy. |
[out] | nfound | The total number of eigenvalues found. 0 <= nfound <= n.
|
[out] | W | The vector W of length n. The first nfound elements contain the selected eigenvalues in ascending order. |
[out] | Z | The n-by-nfound matrix Z, stored in an ldz-by-zcol array.
|
[in] | ldz | The leading dimension of the array Z. ldz >= 1, and if jobz = Vec, ldz >= max(1,n). |
[out] | isuppz | The vector isuppz of length 2*max(1,nfound). The support of the eigenvectors in Z, i.e., the indices indicating the nonzero elements in Z. The i-th eigenvector is nonzero only in elements isuppz( 2*i-1 ) through isuppz( 2*i ). This is an output of lapack::stemr (tridiagonal matrix). The support of the eigenvectors of A is typically 1:n because of the unitary transformations applied by lapack::unmtr . Implemented only for range = All or Index and iu - il = n - 1 |
int64_t lapack::heevr_2stage | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
std::complex< double > * | Z, | ||
int64_t | ldz, | ||
int64_t * | isuppz | ||
) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal.
Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
heevr_2stage
first reduces the matrix A to tridiagonal form T with a call to lapack::hetrd
. Then, whenever possible, heevr_2stage
calls lapack::stemr
to compute eigenspectrum using Relatively Robust Representations. lapack::stemr
computes eigenvalues by the dqds algorithm, while orthogonal eigenvectors are computed from various "good" \(L D L^T\) representations (also known as Relatively Robust Representations). Gram-Schmidt orthogonalization is avoided as far as possible. More specifically, the various steps of the algorithm are as follows.
For each unreduced block (submatrix) of T,
(a) Compute \(T - \sigma I = L D L^T\), so that L and D define all the wanted eigenvalues to high relative accuracy. This means that small relative changes in the entries of D and L cause only small relative changes in the eigenvalues and eigenvectors. The standard (unfactored) representation of the tridiagonal matrix T does not have this property in general.
(b) Compute the eigenvalues to suitable accuracy. If the eigenvectors are desired, the algorithm attains full accuracy of the computed eigenvalues only right before the corresponding vectors have to be computed, see steps c) and d).
(c) For each cluster of close eigenvalues, select a new shift close to the cluster, find a new factorization, and refine the shifted eigenvalues to suitable accuracy.
(d) For each eigenvalue with a large enough relative separation compute the corresponding eigenvector by forming a rank revealing twisted factorization. Go back to (c) for any clusters that remain.
The desired accuracy of the output can be specified by the input parameter abstol.
For more details, see lapack::stemr
documentation and:
Note 1 : heevr_2stage
calls lapack::stemr
when the full spectrum is requested on machines which conform to the ieee-754 floating point standard. heevr_2stage
calls lapack::stebz
and lapack::stein
on non-IEEE machines and when partial spectrum requests are made.
Normal execution of lapack::stemr
may create NaNs and infinities and hence may abort due to a floating point exception in environments which do not handle NaNs and infinities in the IEEE standard default manner.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::syevr_2stage
.
[in] | jobz |
|
[in] | range |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | vl | If range=Value, the lower bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | vu | If range=Value, the upper bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | il | If range=Index, the index of the smallest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | iu | If range=Index, the index of the largest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | abstol | The absolute error tolerance for the eigenvalues. 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 + eps * max( |a|,|b| ), where eps is the machine precision. If abstol is less than or equal to zero, then eps*|T| will be used in its place, where |T| is the 1-norm of the tridiagonal matrix obtained by reducing A to tridiagonal form. See "Computing Small Singular Values of Bidiagonal Matrices with Guaranteed High Relative Accuracy," by Demmel and Kahan, LAPACK Working Note #3. If high relative accuracy is important, set abstol to DLAMCH( 'Safe minimum' ). Doing so will guarantee that eigenvalues are computed to high relative accuracy when possible in future releases. The current code does not make any guarantees about high relative accuracy, but future releases will. See J. Barlow and J. Demmel, "Computing Accurate Eigensystems of Scaled Diagonally Dominant Matrices", LAPACK Working Note #7, for a discussion of which matrices define their eigenvalues to high relative accuracy. |
[out] | nfound | The total number of eigenvalues found. 0 <= nfound <= n.
|
[out] | W | The vector W of length n. The first nfound elements contain the selected eigenvalues in ascending order. |
[out] | Z | The n-by-nfound matrix Z, stored in an ldz-by-zcol array.
|
[in] | ldz | The leading dimension of the array Z. ldz >= 1, and if jobz = Vec, ldz >= max(1,n). |
[out] | isuppz | The vector isuppz of length 2*max(1,nfound). The support of the eigenvectors in Z, i.e., the indices indicating the nonzero elements in Z. The i-th eigenvector is nonzero only in elements isuppz( 2*i-1 ) through isuppz( 2*i ). This is an output of lapack::stemr (tridiagonal matrix). The support of the eigenvectors of A is typically 1:n because of the unitary transformations applied by lapack::unmtr . Implemented only for range = All or Index and iu - il = n - 1 |
All details about the 2-stage techniques are available in:
Azzam Haidar, Hatem Ltaief, and Jack Dongarra. Parallel reduction to condensed forms for symmetric eigenvalue problems using aggregated fine-grained and memory-aware kernels. In Proceedings of 2011 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '11), New York, NY, USA, Article 8, 11 pages. http://doi.acm.org/10.1145/2063384.2063394
A. Haidar, J. Kurzak, P. Luszczek, 2013. An improved parallel singular value algorithm and its implementation for multicore hardware, In Proceedings of 2013 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '13). Denver, Colorado, USA, 2013. Article 90, 12 pages. http://doi.acm.org/10.1145/2503210.2503292
A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. A novel hybrid CPU-GPU generalized eigensolver for electronic structure calculations based on fine-grained memory aware tasks. International Journal of High Performance Computing Applications. Volume 28 Issue 2, Pages 196-209, May 2014. http://hpc.sagepub.com/content/28/2/196
int64_t lapack::heevx | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
std::complex< double > * | Z, | ||
int64_t | ldz, | ||
int64_t * | ifail | ||
) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A.
Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | jobz |
|
[in] | range |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | vl | If range=Value, the lower bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | vu | If range=Value, the upper bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | il | If range=Index, the index of the smallest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | iu | If range=Index, the index of the largest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | abstol | The absolute error tolerance for the eigenvalues. 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 + eps * max(|a|, |b|), where eps is the machine precision. If abstol is less than or equal to zero, then eps*|T| will be used in its place, where |T| is the 1-norm of 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*DLAMCH('S'), not zero. If this routine returns with return value > 0, indicating that some eigenvectors did not converge, try setting abstol to 2*DLAMCH('S'). See "Computing Small Singular Values of Bidiagonal Matrices with Guaranteed High Relative Accuracy," by Demmel and Kahan, LAPACK Working Note #3. |
[out] | nfound | The total number of eigenvalues found. 0 <= nfound <= n.
|
[out] | W | The vector W of length n. On normal exit, the first nfound elements contain the selected eigenvalues in ascending order. |
[out] | Z | The n-by-nfound matrix Z, stored in an ldz-by-zcol array.
|
[in] | ldz | The leading dimension of the array Z. ldz >= 1, and if jobz = Vec, ldz >= max(1,n). |
[out] | ifail | The vector ifail of length n.
|
int64_t lapack::heevx_2stage | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
std::complex< double > * | Z, | ||
int64_t | ldz, | ||
int64_t * | ifail | ||
) |
Computes selected eigenvalues and, optionally, eigenvectors of a Hermitian matrix A using the 2-stage technique for the reduction to tridiagonal.
Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | jobz |
|
[in] | range |
|
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | vl | If range=Value, the lower bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | vu | If range=Value, the upper bound of the interval to be searched for eigenvalues. vl < vu. Not referenced if range = All or Index. |
[in] | il | If range=Index, the index of the smallest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | iu | If range=Index, the index of the largest eigenvalue to be returned. 1 <= il <= iu <= n, if n > 0; il = 1 and iu = 0 if n = 0. Not referenced if range = All or Value. |
[in] | abstol | The absolute error tolerance for the eigenvalues. 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 + eps * max(|a|, |b|), where eps is the machine precision. If abstol is less than or equal to zero, then eps*|T| will be used in its place, where |T| is the 1-norm of 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*DLAMCH('S'), not zero. If this routine returns with return value > 0, indicating that some eigenvectors did not converge, try setting abstol to 2*DLAMCH('S'). See "Computing Small Singular Values of Bidiagonal Matrices with Guaranteed High Relative Accuracy," by Demmel and Kahan, LAPACK Working Note #3. |
[out] | nfound | The total number of eigenvalues found. 0 <= nfound <= n.
|
[out] | W | The vector W of length n. On normal exit, the first nfound elements contain the selected eigenvalues in ascending order. |
[out] | Z | The n-by-nfound matrix Z, stored in an ldz-by-zcol array.
|
[in] | ldz | The leading dimension of the array Z. ldz >= 1, and if jobz = Vec, ldz >= max(1,n). |
[out] | ifail | The vector ifail of length n.
|
All details about the 2-stage techniques are available in:
Azzam Haidar, Hatem Ltaief, and Jack Dongarra. Parallel reduction to condensed forms for symmetric eigenvalue problems using aggregated fine-grained and memory-aware kernels. In Proceedings of 2011 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '11), New York, NY, USA, Article 8, 11 pages. http://doi.acm.org/10.1145/2063384.2063394
A. Haidar, J. Kurzak, P. Luszczek, 2013. An improved parallel singular value algorithm and its implementation for multicore hardware, In Proceedings of 2013 International Conference for High Performance Computing, Networking, Storage and Analysis (SC '13). Denver, Colorado, USA, 2013. Article 90, 12 pages. http://doi.acm.org/10.1145/2503210.2503292
A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. A novel hybrid CPU-GPU generalized eigensolver for electronic structure calculations based on fine-grained memory aware tasks. International Journal of High Performance Computing Applications. Volume 28 Issue 2, Pages 196-209, May 2014. http://hpc.sagepub.com/content/28/2/196
int64_t lapack::syev | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
int64_t lapack::syev_2stage | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
int64_t lapack::syevd | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
int64_t lapack::syevd_2stage | ( | lapack::Job | jobz, |
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double * | W | ||
) |
int64_t lapack::syevr | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
double * | Z, | ||
int64_t | ldz, | ||
int64_t * | isuppz | ||
) |
int64_t lapack::syevr_2stage | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
double * | Z, | ||
int64_t | ldz, | ||
int64_t * | isuppz | ||
) |
int64_t lapack::syevx | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
double * | Z, | ||
int64_t | ldz, | ||
int64_t * | ifail | ||
) |
int64_t lapack::syevx_2stage | ( | lapack::Job | jobz, |
lapack::Range | range, | ||
lapack::Uplo | uplo, | ||
int64_t | n, | ||
double * | A, | ||
int64_t | lda, | ||
double | vl, | ||
double | vu, | ||
int64_t | il, | ||
int64_t | iu, | ||
double | abstol, | ||
int64_t * | nfound, | ||
double * | W, | ||
double * | Z, | ||
int64_t | ldz, | ||
int64_t * | ifail | ||
) |