Spectra 1.1.0
Header-only C++ Library for Large Scale Eigenvalue Problems
Loading...
Searching...
No Matches
Spectra::HermEigsSolver< OpType > Class Template Reference

#include <Spectra/HermEigsSolver.h>

Inheritance diagram for Spectra::HermEigsSolver< OpType >:
Spectra::HermEigsBase< OpType, BOpType >

Public Member Functions

 HermEigsSolver (OpType &op, Index nev, Index ncv)
 
- Public Member Functions inherited from Spectra::HermEigsBase< OpType, BOpType >
void init (const Scalar *init_resid)
 
void init ()
 
Index compute (SortRule selection=SortRule::LargestMagn, Index maxit=1000, RealScalar tol=1e-10, SortRule sorting=SortRule::LargestAlge)
 
CompInfo info () const
 
Index num_iterations () const
 
Index num_operations () const
 
RealVector eigenvalues () const
 
virtual Matrix eigenvectors (Index nvec) const
 
virtual Matrix eigenvectors () const
 

Detailed Description

template<typename OpType = DenseHermMatProd<double>>
class Spectra::HermEigsSolver< OpType >

This class implements the eigen solver for Hermitian matrices, i.e., to solve \(Ax=\lambda x\) where \(A\) 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
OpTypeThe 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 <Spectra/HermEigsSolver.h>
// <Spectra/MatOp/DenseHermMatProd.h> is implicitly included
#include <iostream>
using namespace Spectra;
int main()
{
// We are going to calculate the eigenvalues of M
Eigen::MatrixXcd A = Eigen::MatrixXcd::Random(10, 10);
Eigen::MatrixXcd M = A + A.adjoint();
// Construct matrix operation object using the wrapper class DenseHermMatProd
OpType op(M);
// Construct eigen solver object, requesting the largest three eigenvalues
HermEigsSolver<OpType> eigs(op, 3, 6);
// Initialize and compute
eigs.init();
int nconv = eigs.compute(SortRule::LargestAlge);
// Retrieve results
// Eigenvalues are real-valued, and eigenvectors are complex-valued
Eigen::VectorXd evalues;
if (eigs.info() == CompInfo::Successful)
evalues = eigs.eigenvalues();
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
Eigen::MatrixXcd evecs = eigs.eigenvectors();
std::cout << "Eigenvectors:\n" << evecs << std::endl;
return 0;
}
HermEigsSolver(OpType &op, Index nev, Index ncv)
@ Successful
Computation was successful.
Definition CompInfo.h:19

And here is an example for user-supplied matrix operation class.

#include <Eigen/Core>
#include <Spectra/HermEigsSolver.h>
#include <iostream>
using namespace Spectra;
// M = diag(1+0i, 2+0i, ..., 10+0i)
class MyDiagonalTen
{
public:
using Scalar = std::complex<double>; // A typedef named "Scalar" is required
int rows() const { return 10; }
int cols() const { return 10; }
// y_out = M * x_in
void perform_op(const Scalar *x_in, Scalar *y_out) const
{
for (int i = 0; i < rows(); i++)
{
y_out[i] = x_in[i] * Scalar(i + 1, 0);
}
}
};
int main()
{
MyDiagonalTen op;
eigs.init();
eigs.compute(SortRule::LargestAlge);
if (eigs.info() == CompInfo::Successful)
{
Eigen::VectorXd evalues = eigs.eigenvalues();
// Will get (10, 9, 8)
std::cout << "Eigenvalues found:\n" << evalues << std::endl;
Eigen::MatrixXcd evecs = eigs.eigenvectors();
std::cout << "Eigenvectors:\n" << evecs << std::endl;
}
return 0;
}

Definition at line 122 of file HermEigsSolver.h.

Constructor & Destructor Documentation

◆ HermEigsSolver()

template<typename OpType = DenseHermMatProd<double>>
Spectra::HermEigsSolver< OpType >::HermEigsSolver ( OpType & op,
Index nev,
Index ncv )
inline

Constructor to create a solver object.

Parameters
opThe matrix operation object that implements the matrix-vector multiplication operation of \(A\): calculating \(Av\) for any vector \(v\). Users could either create the object from the wrapper class such as DenseHermMatProd, or define their own that implements all the public members as in DenseHermMatProd.
nevNumber of eigenvalues requested. This should satisfy \(1\le nev \le n-1\), where \(n\) is the size of matrix.
ncvParameter that controls the convergence speed of the algorithm. Typically a larger ncv means faster convergence, but it may also result in greater memory use and more matrix operations in each iteration. This parameter must satisfy \(nev < ncv \le n\), and is advised to take \(ncv \ge 2\cdot nev\).

Definition at line 145 of file HermEigsSolver.h.


The documentation for this class was generated from the following file: