MFEM v4.7.0
Finite element discretization library
|
#include <slepc.hpp>
Public Types | |
enum | Which { LARGEST_MAGNITUDE , SMALLEST_MAGNITUDE , LARGEST_REAL , SMALLEST_REAL , LARGEST_IMAGINARY , SMALLEST_IMAGINARY , TARGET_MAGNITUDE , TARGET_REAL } |
Target spectrum for the eigensolver. More... | |
enum | SpectralTransformation { SHIFT , SHIFT_INVERT } |
Spectral transformations that can be used by the solver in order to accelerate the convergence to the target eignevalues. More... | |
Public Member Functions | |
SlepcEigenSolver (MPI_Comm comm, const std::string &prefix=std::string()) | |
Constructors. | |
virtual | ~SlepcEigenSolver () |
void | SetTol (real_t tol) |
Set solver convergence tolerance relative to the magnitude of the eigenvalue. | |
void | SetMaxIter (int max_iter) |
Set maximum number of iterations allowed in the call to SlepcEigenSolver::Solve. | |
void | SetNumModes (int num_eigs) |
Set the number of eigenmodes to compute. | |
void | SetOperator (const PetscParMatrix &op) |
Set operator for standard eigenvalue problem. | |
void | SetOperators (const PetscParMatrix &op, const PetscParMatrix &opB) |
Set operators for generalized eigenvalue problem. | |
void | Customize (bool customize=true) const |
Customize object with options set. | |
void | Solve () |
Solve the eigenvalue problem for the specified number of eigenvalues. | |
int | GetNumConverged () |
Get the number of converged eigenvalues after the call to SlepcEigenSolver::Solve. | |
void | GetEigenvalue (unsigned int i, real_t &lr) const |
Get the ith eigenvalue after the system has been solved. | |
void | GetEigenvalue (unsigned int i, real_t &lr, real_t &lc) const |
Get the ith eigenvalue after the system has been solved. | |
void | GetEigenvector (unsigned int i, Vector &vr) const |
Get the ith eigenvector after the system has been solved. | |
void | GetEigenvector (unsigned int i, Vector &vr, Vector &vc) const |
Get the ith eigenvector after the system has been solved. | |
void | SetWhichEigenpairs (Which which) |
Set the which eigenvalues the solver will target and the order they will be indexed in. | |
void | SetTarget (real_t target) |
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or SlepcEigenSolver::TARGET_REAL in the SlepcEigenSolver::SetWhichEigenpairs method. | |
void | SetSpectralTransformation (SpectralTransformation transformation) |
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT and SlepcEigenSolver::SHIFT_INVERT are available. | |
operator slepc::EPS () const | |
Conversion function to SLEPc's EPS type. | |
operator PetscObject () const | |
Conversion function to PetscObject. | |
Target spectrum for the eigensolver.
This will define the order in which the eigenvalues/eigenvectors are indexed after the call to SlepcEigenSolver::Solve.
mfem::SlepcEigenSolver::SlepcEigenSolver | ( | MPI_Comm | comm, |
const std::string & | prefix = std::string() ) |
void mfem::SlepcEigenSolver::Customize | ( | bool | customize = true | ) | const |
void mfem::SlepcEigenSolver::GetEigenvalue | ( | unsigned int | i, |
real_t & | lr ) const |
Get the ith eigenvalue after the system has been solved.
[in] | i | The index for the eigenvalue you want ordered by SlepcEigenSolver::SetWhichEigenpairs |
[out] | lr | The real component of the eigenvalue |
Get the ith eigenvalue after the system has been solved.
[in] | i | The index for the eigenvalue you want ordered by SlepcEigenSolver::SetWhichEigenpairs |
[out] | lr | The real component of the eigenvalue |
[out] | lc | The imaginary component of the eigenvalue |
void mfem::SlepcEigenSolver::GetEigenvector | ( | unsigned int | i, |
Vector & | vr ) const |
Get the ith eigenvector after the system has been solved.
[in] | i | The index for the eigenvector you want ordered by SlepcEigenSolver::SetWhichEigenpairs |
[out] | vr | The real components of the eigenvector |
Get the ith eigenvector after the system has been solved.
[in] | i | The index for the eigenvector you want ordered by SlepcEigenSolver::SetWhichEigenpairs |
[out] | vr | The real components of the eigenvector |
[out] | vc | The imaginary components of the eigenvector |
int mfem::SlepcEigenSolver::GetNumConverged | ( | ) |
Get the number of converged eigenvalues after the call to SlepcEigenSolver::Solve.
|
inline |
|
inline |
void mfem::SlepcEigenSolver::SetMaxIter | ( | int | max_iter | ) |
Set maximum number of iterations allowed in the call to SlepcEigenSolver::Solve.
void mfem::SlepcEigenSolver::SetNumModes | ( | int | num_eigs | ) |
void mfem::SlepcEigenSolver::SetOperator | ( | const PetscParMatrix & | op | ) |
void mfem::SlepcEigenSolver::SetOperators | ( | const PetscParMatrix & | op, |
const PetscParMatrix & | opB ) |
void mfem::SlepcEigenSolver::SetSpectralTransformation | ( | SlepcEigenSolver::SpectralTransformation | transformation | ) |
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT and SlepcEigenSolver::SHIFT_INVERT are available.
void mfem::SlepcEigenSolver::SetTarget | ( | real_t | target | ) |
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or SlepcEigenSolver::TARGET_REAL in the SlepcEigenSolver::SetWhichEigenpairs method.
void mfem::SlepcEigenSolver::SetTol | ( | real_t | tol | ) |
void mfem::SlepcEigenSolver::SetWhichEigenpairs | ( | SlepcEigenSolver::Which | which | ) |
Set the which eigenvalues the solver will target and the order they will be indexed in.
For SlepcEigenSolver::TARGET_MAGNITUDE or SlepcEigenSolver::TARGET_REAL you will also need to set the target value with SlepcEigenSolver::SetTarget.
void mfem::SlepcEigenSolver::Solve | ( | ) |