DOLFIN
DOLFIN C++ interface
|
#include <SLEPcEigenSolver.h>
Public Member Functions | |
SLEPcEigenSolver (MPI_Comm comm) | |
Create eigenvalue solver. | |
SLEPcEigenSolver (EPS eps) | |
Create eigenvalue solver from EPS object. | |
SLEPcEigenSolver (std::shared_ptr< const PETScMatrix > A) | |
Create eigenvalue solver for Ax = . | |
SLEPcEigenSolver (MPI_Comm comm, std::shared_ptr< const PETScMatrix > A) | |
Create eigenvalue solver for Ax = x. | |
SLEPcEigenSolver (std::shared_ptr< const PETScMatrix > A, std::shared_ptr< const PETScMatrix > B) | |
Create eigenvalue solver for Ax = x on MPI_COMM_WORLD. | |
SLEPcEigenSolver (MPI_Comm comm, std::shared_ptr< const PETScMatrix > A, std::shared_ptr< const PETScMatrix > B) | |
Create eigenvalue solver for Ax = x. | |
~SLEPcEigenSolver () | |
Destructor. | |
void | set_operators (std::shared_ptr< const PETScMatrix > A, std::shared_ptr< const PETScMatrix > B) |
void | solve () |
Compute all eigenpairs of the matrix A (solve Ax = x) | |
void | solve (std::size_t n) |
Compute the n first eigenpairs of the matrix A (solve Ax = x) | |
void | get_eigenvalue (double &lr, double &lc, std::size_t i) const |
Get ith eigenvalue. | |
void | get_eigenpair (double &lr, double &lc, GenericVector &r, GenericVector &c, std::size_t i) const |
Get ith eigenpair. | |
void | get_eigenpair (double &lr, double &lc, PETScVector &r, PETScVector &c, std::size_t i) const |
Get ith eigenpair. | |
std::size_t | get_iteration_number () const |
Get the number of iterations used by the solver. | |
std::size_t | get_number_converged () const |
Get the number of converged eigenvalues. | |
void | set_deflation_space (const VectorSpaceBasis &deflation_space) |
void | set_initial_space (const VectorSpaceBasis &initial_space) |
void | set_options_prefix (std::string options_prefix) |
std::string | get_options_prefix () const |
void | set_from_options () const |
Set options from PETSc options database. | |
EPS | eps () const |
Return SLEPc EPS pointer. | |
Public Member Functions inherited from dolfin::Variable | |
Variable () | |
Create unnamed variable. | |
Variable (const std::string name, const std::string label) | |
Create variable with given name and label. | |
Variable (const Variable &variable) | |
Copy constructor. | |
virtual | ~Variable () |
Destructor. | |
const Variable & | operator= (const Variable &variable) |
Assignment operator. | |
void | rename (const std::string name, const std::string label) |
Rename variable. | |
std::string | name () const |
Return name. | |
std::string | label () const |
Return label (description) | |
std::size_t | id () const |
virtual std::string | str (bool verbose) const |
Return informal string representation (pretty-print) | |
Public Member Functions inherited from dolfin::PETScObject | |
PETScObject () | |
Constructor. Ensures that PETSc has been initialised. | |
virtual | ~PETScObject () |
Destructor. | |
Static Public Member Functions | |
static Parameters | default_parameters () |
Default parameter values. | |
Static Public Member Functions inherited from dolfin::PETScObject | |
static void | petsc_error (int error_code, std::string filename, std::string petsc_function) |
Print error message for PETSc calls that return an error. | |
Additional Inherited Members | |
Public Attributes inherited from dolfin::Variable | |
Parameters | parameters |
Parameters. | |
This class provides an eigenvalue solver for PETSc matrices. It is a wrapper for the SLEPc eigenvalue solver.
The following parameters may be specified to control the solver.
This parameter controls which part of the spectrum to compute. Possible values are
"largest magnitude" (eigenvalues with largest magnitude) "smallest magnitude" (eigenvalues with smallest magnitude) "largest real" (eigenvalues with largest double part) "smallest real" (eigenvalues with smallest double part) "largest imaginary" (eigenvalues with largest imaginary part) "smallest imaginary" (eigenvalues with smallest imaginary part) "target magnitude" (eigenvalues closest to target in magnitude) "target real" (eigenvalues closest to target in real part) "target imaginary" (eigenvalues closest to target in imaginary part)
This parameter controls which algorithm is used by SLEPc. Possible values are
"power" (power iteration) "subspace" (subspace iteration) "arnoldi" (Arnoldi) "lanczos" (Lanczos) "krylov-schur" (Krylov-Schur) "lapack" (LAPACK, all values, direct, small systems only) "arpack" (ARPACK)
This parameter controls the tolerance used by SLEPc. Possible values are positive double numbers.
This parameter controls the maximum number of iterations used by SLEPc. Possible values are positive integers.
Note that both the tolerance and the number of iterations must be specified if either one is specified.
This parameter can be used to give extra information about the type of the eigenvalue problem. Some solver types require this extra piece of information. Possible values are:
"hermitian" (Hermitian) "non_hermitian" (Non-Hermitian) "gen_hermitian" (Generalized Hermitian) "gen_non_hermitian" (Generalized Non-Hermitian) "pos_gen_non_hermitian" (Generalized Non-Hermitian with positive semidefinite B)
This parameter controls the application of a spectral transform. A spectral transform can be used to enhance the convergence of the eigensolver and in particular to only compute eigenvalues in the interior of the spectrum. Possible values are:
"shift-and-invert" (A shift-and-invert transform)
Note that if a spectral transform is given, then also a non-zero spectral shift parameter has to be provided.
The default is no spectral transform.
This parameter controls the spectral shift used by the spectral transform and must be provided if a spectral transform is given. The possible values are real numbers.
std::string SLEPcEigenSolver::get_options_prefix | ( | ) | const |
Returns the prefix used by PETSc when searching the PETSc options database
void SLEPcEigenSolver::set_deflation_space | ( | const VectorSpaceBasis & | deflation_space | ) |
Set deflation space. The VectorSpaceBasis does not need to be orthonormal.
void SLEPcEigenSolver::set_initial_space | ( | const VectorSpaceBasis & | initial_space | ) |
Set inital space. The VectorSpaceBasis does not need to be orthonormal.
void SLEPcEigenSolver::set_operators | ( | std::shared_ptr< const PETScMatrix > | A, |
std::shared_ptr< const PETScMatrix > | B | ||
) |
Set opeartors (B may be nullptr for regular eigenvalues problems)
void SLEPcEigenSolver::set_options_prefix | ( | std::string | options_prefix | ) |
Sets the prefix used by PETSc when searching the PETSc options database