Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/slepc_solver.h>
Public Member Functions | |
SolverBase (SolverControl &cn, const MPI_Comm &mpi_communicator) | |
virtual | ~SolverBase () |
template<typename OutputVector > | |
void | solve (const PETScWrappers::MatrixBase &A, std::vector< PetscScalar > &eigenvalues, std::vector< OutputVector > &eigenvectors, const unsigned int n_eigenpairs=1) |
template<typename OutputVector > | |
void | solve (const PETScWrappers::MatrixBase &A, const PETScWrappers::MatrixBase &B, std::vector< PetscScalar > &eigenvalues, std::vector< OutputVector > &eigenvectors, const unsigned int n_eigenpairs=1) |
template<typename OutputVector > | |
void | solve (const PETScWrappers::MatrixBase &A, const PETScWrappers::MatrixBase &B, std::vector< double > &real_eigenvalues, std::vector< double > &imag_eigenvalues, std::vector< OutputVector > &real_eigenvectors, std::vector< OutputVector > &imag_eigenvectors, const unsigned int n_eigenpairs=1) |
void | set_initial_vector (const PETScWrappers::VectorBase &this_initial_vector) 1 |
template<typename Vector > | |
void | set_initial_space (const std::vector< Vector > &initial_space) |
void | set_transformation (SLEPcWrappers::TransformationBase &this_transformation) |
void | set_target_eigenvalue (const PetscScalar &this_target) |
void | set_which_eigenpairs (EPSWhich set_which) |
void | set_problem_type (EPSProblemType set_problem) |
void | get_solver_state (const SolverControl::State state) |
SolverControl & | control () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcSLEPcWrappersUsageError () |
static ::ExceptionBase & | ExcSLEPcError (int arg1) |
static ::ExceptionBase & | ExcSLEPcEigenvectorConvergenceMismatchError (int arg1, int arg2) |
Protected Member Functions | |
void | solve (const unsigned int n_eigenpairs, unsigned int *n_converged) |
void | get_eigenpair (const unsigned int index, PetscScalar &eigenvalues, PETScWrappers::VectorBase &eigenvectors) |
void | get_eigenpair (const unsigned int index, double &real_eigenvalues, double &imag_eigenvalues, PETScWrappers::VectorBase &real_eigenvectors, PETScWrappers::VectorBase &imag_eigenvectors) |
void | set_matrices (const PETScWrappers::MatrixBase &A) |
void | set_matrices (const PETScWrappers::MatrixBase &A, const PETScWrappers::MatrixBase &B) |
Protected Attributes | |
SolverControl & | solver_control |
const MPI_Comm | mpi_communicator |
EPS | eps |
Static Private Member Functions | |
static int | convergence_test (EPS eps, PetscScalar real_eigenvalue, PetscScalar imag_eigenvalue, PetscReal residual_norm, PetscReal *estimated_error, void *solver_control) |
Private Attributes | |
EPSConvergedReason | reason |
Base class for solver classes using the SLEPc solvers. Since solvers in SLEPc are selected based on flags passed to a generic solver object, basically all the actual solver calls happen in this class, and derived classes simply set the right flags to select one solver or another, or to set certain parameters for individual solvers.
Definition at line 140 of file slepc_solver.h.
SLEPcWrappers::SolverBase::SolverBase | ( | SolverControl & | cn, |
const MPI_Comm & | mpi_communicator | ||
) |
Constructor. Takes the MPI communicator over which parallel computations are to happen.
Definition at line 36 of file slepc_solver.cc.
|
virtual |
Destructor.
Definition at line 60 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::solve | ( | const PETScWrappers::MatrixBase & | A, |
std::vector< PetscScalar > & | eigenvalues, | ||
std::vector< OutputVector > & | eigenvectors, | ||
const unsigned int | n_eigenpairs = 1 |
||
) |
Composite method that solves the eigensystem \(Ax=\lambda x\). The eigenvector sent in has to have at least one element that we can use as a template when resizing, since we do not know the parameters of the specific vector class used (i.e. local_dofs for MPI vectors). However, while copying eigenvectors, at least twice the memory size of eigenvectors
is being used (and can be more). To avoid doing this, the fairly standard calling sequence executed here is used: Set up matrices for solving; Actually solve the system; Gather the solution(s).
This is declared here to make it possible to take a std::vector of different PETScWrappers vector types
Definition at line 672 of file slepc_solver.h.
void SLEPcWrappers::SolverBase::solve | ( | const PETScWrappers::MatrixBase & | A, |
const PETScWrappers::MatrixBase & | B, | ||
std::vector< PetscScalar > & | eigenvalues, | ||
std::vector< OutputVector > & | eigenvectors, | ||
const unsigned int | n_eigenpairs = 1 |
||
) |
Same as above, but here a composite method for solving the system \(A x=\lambda B x\), for real matrices, vectors, and values \(A, B, x, \lambda\).
Definition at line 703 of file slepc_solver.h.
void SLEPcWrappers::SolverBase::solve | ( | const PETScWrappers::MatrixBase & | A, |
const PETScWrappers::MatrixBase & | B, | ||
std::vector< double > & | real_eigenvalues, | ||
std::vector< double > & | imag_eigenvalues, | ||
std::vector< OutputVector > & | real_eigenvectors, | ||
std::vector< OutputVector > & | imag_eigenvectors, | ||
const unsigned int | n_eigenpairs = 1 |
||
) |
Same as above, but here a composite method for solving the system \(A x=\lambda B x\) with real matrices \(A, B\) and imaginary eigenpairs \(x, \lambda\).
Definition at line 740 of file slepc_solver.h.
void SLEPcWrappers::SolverBase::set_initial_vector | ( | const PETScWrappers::VectorBase & | this_initial_vector | ) |
Set the initial vector for the solver.
Definition at line 101 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::set_initial_space | ( | const std::vector< Vector > & | initial_space | ) |
Set the initial vector space for the solver.
By default, SLEPc initializes the starting vector or the initial subspace randomly.
Definition at line 790 of file slepc_solver.h.
void SLEPcWrappers::SolverBase::set_transformation | ( | SLEPcWrappers::TransformationBase & | this_transformation | ) |
Set the spectral transformation to be used.
Definition at line 92 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::set_target_eigenvalue | ( | const PetscScalar & | this_target | ) |
Set target eigenvalues in the spectrum to be computed. By default, no target is set.
Definition at line 116 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::set_which_eigenpairs | ( | EPSWhich | set_which | ) |
Indicate which part of the spectrum is to be computed. By default largest magnitude eigenvalues are computed.
Definition at line 126 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::set_problem_type | ( | EPSProblemType | set_problem | ) |
Specify the type of the eigenspectrum problem. This can be used to exploit known symmetries of the matrices that make up the standard/generalized eigenspectrum problem. By default a non-Hermitian problem is assumed.
Definition at line 134 of file slepc_solver.cc.
void SLEPcWrappers::SolverBase::get_solver_state | ( | const SolverControl::State | state | ) |
Take the information provided from SLEPc and checks it against deal.II's own SolverControl objects to see if convergence has been reached.
Definition at line 270 of file slepc_solver.cc.
SolverControl & SLEPcWrappers::SolverBase::control | ( | ) | const |
Access to the object that controls convergence.
Definition at line 296 of file slepc_solver.cc.
|
protected |
Solve the linear system for n_eigenpairs
eigenstates. Parameter n_converged
contains the actual number of eigenstates that have converged; this can be both fewer or more than n_eigenpairs, depending on the SLEPc eigensolver used.
Definition at line 141 of file slepc_solver.cc.
|
protected |
Access the real parts of solutions for a solved eigenvector problem, pair index solutions, \(\text{index}\,\in\,0\hdots \text{n\_converged}-1\).
Definition at line 228 of file slepc_solver.cc.
|
protected |
Access the real and imaginary parts of solutions for a solved eigenvector problem, pair index solutions, \(\text{index}\,\in\,0\hdots \text{n\_converged}-1\).
Definition at line 241 of file slepc_solver.cc.
|
protected |
Initialize solver for the linear system \(Ax=\lambda x\). (Note: this is required before calling solve ())
Definition at line 75 of file slepc_solver.cc.
|
protected |
Same as above, but here initialize solver for the linear system \(A x=\lambda B x\).
Definition at line 83 of file slepc_solver.cc.
|
staticprivate |
A function that can be used in SLEPc as a callback to check on convergence.
Definition at line 302 of file slepc_solver.cc.
|
protected |
Reference to the object that controls convergence of the iterative solver.
Definition at line 299 of file slepc_solver.h.
|
protected |
Copy of the MPI communicator object to be used for the solver.
Definition at line 304 of file slepc_solver.h.
|
protected |
Objects for Eigenvalue Problem Solver.
Definition at line 358 of file slepc_solver.h.
|
private |
Convergence reason.
Definition at line 364 of file slepc_solver.h.