SLATE 2024.05.31
Software for Linear Algebra Targeting Exascale
|
\(Ax = \lambda x\) More...
Functions | |
template<typename scalar_t > | |
void | slate::heev (HermitianMatrix< scalar_t > &A, std::vector< blas::real_type< scalar_t > > &Lambda, Matrix< scalar_t > &Z, Options const &opts) |
Distributed parallel Hermitian matrix eigen decomposition. | |
\(Ax = \lambda x\)
void slate::heev | ( | HermitianMatrix< scalar_t > & | A, |
std::vector< blas::real_type< scalar_t > > & | Lambda, | ||
Matrix< scalar_t > & | Z, | ||
Options const & | opts | ||
) |
Distributed parallel Hermitian matrix eigen decomposition.
heev Computes all eigenvalues and, optionally, eigenvectors of a Hermitian matrix A. The matrix A is preliminary reduced to tridiagonal form using a two-stage approach: First stage: reduction to band tridiagonal form (see he2hb); Second stage: reduction from band to tridiagonal form (see hb2st).
scalar_t | One of float, double, std::complex<float>, std::complex<double>. |
[in] | A | On entry, the n-by-n Hermitian matrix \(A\). On exit, contents are destroyed. |
[out] | Lambda | The vector Lambda of length n. If successful, the eigenvalues in ascending order. |
[out] | Z | On entry, if Z is empty, does not compute eigenvectors. Otherwise, the n-by-n matrix \(Z\) to store eigenvectors. On exit, orthonormal eigenvectors of the matrix A. |
[in] | opts | Additional options, as map of name = value pairs. Possible options:
|