template<typename OpType = DenseHermMatProd<double>>
class Spectra::HermEigsSolver< OpType >
This class implements the eigen solver for Hermitian matrices, i.e., to solve
where
is Hermitian. An Hermitian matrix is a complex square matrix that is equal to its own conjugate transpose. It is known that all Hermitian matrices have real-valued eigenvalues.
- Template Parameters
-
OpType | The name of the matrix operation class. Users could either use the wrapper classes such as DenseHermMatProd and SparseHermMatProd, or define their own that implements the type definition Scalar and all the public member functions as in DenseHermMatProd. |
Below is an example that demonstrates the usage of this class.
#include <Eigen/Core>
#include <iostream>
{
Eigen::MatrixXcd
A = Eigen::MatrixXcd::Random(10, 10);
Eigen::MatrixXcd
M =
A +
A.adjoint();
eigs.init();
Eigen::VectorXd evalues;
evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
return 0;
}
And here is an example for user-supplied matrix operation class.
#include <Eigen/Core>
#include <iostream>
class MyDiagonalTen
{
public:
using Scalar = std::complex<double>;
int rows()
const {
return 10; }
int cols()
const {
return 10; }
{
{
}
}
};
{
MyDiagonalTen op;
eigs.init();
{
Eigen::VectorXd evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
}
return 0;
}
Definition at line 116 of file HermEigsSolver.h.