SLATE 2024.05.31
Software for Linear Algebra Targeting Exascale
|
Functions | |
template<typename scalar_t > | |
void | slate::hegst (int64_t itype, HermitianMatrix< scalar_t > &A, HermitianMatrix< scalar_t > &B, Options const &opts) |
Distributed parallel reduction of a complex Hermitian positive-definite generalized eigenvalue problem to the standard form. | |
void slate::hegst | ( | int64_t | itype, |
HermitianMatrix< scalar_t > & | A, | ||
HermitianMatrix< scalar_t > & | B, | ||
Options const & | opts | ||
) |
Distributed parallel reduction of a complex Hermitian positive-definite generalized eigenvalue problem to the standard form.
Reduces a complex Hermitian positive-definite generalized eigenvalue problem to standard form, as follows:
itype | Problem |
---|---|
itype = 1 | \(A x = \lambda B x\) |
itype = 2 | \(A B x = \lambda x\) |
itype = 3 | \(B A x = \lambda x\) |
Before calling slate::hegst
, you must call slate::potrf
to compute the Cholesky factorization: \(B = L L^H\) or \(B = U^H U\).
scalar_t | One of float, double, std::complex<float>, std::complex<double>. |
[in] | itype |
|
[in,out] | A | On entry, the n-by-n Hermitian matrix \(A\). On exit, the upper or lower triangle is overwritten by the upper or lower triangle of C, as follows:
|
[in] | B | On entry, the triangular factor from the Cholesky factorization of \(B\), as returned by |slatepotrf|. |
[in] | opts | Additional options, as map of name = value pairs. Possible options:
|
TODO: return value