Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | Private Attributes | List of all members
SLEPcWrappers::SolverBase Class Reference

#include <deal.II/lac/slepc_solver.h>

Inheritance diagram for SLEPcWrappers::SolverBase:
[legend]

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)
 
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)
 
SolverControlcontrol () 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

SolverControlsolver_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
 

Detailed Description

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 146 of file slepc_solver.h.

Constructor & Destructor Documentation

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 35 of file slepc_solver.cc.

SLEPcWrappers::SolverBase::~SolverBase ( )
virtual

Destructor.

Definition at line 60 of file slepc_solver.cc.

Member Function Documentation

template<typename OutputVector >
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).

Note
Note that the number of converged eigenvectors can be larger than the number of eigenvectors requested; this is due to a round off error (success) of the eigenproblem solver context. If this is found to be the case we simply do not bother with more eigenpairs than requested, but handle that it may be more than specified by ignoring any extras. By default one eigenvector/eigenvalue pair is computed.

This is declared here to make it possible to take a std::vector of different PETScWrappers vector types

Definition at line 668 of file slepc_solver.h.

template<typename OutputVector >
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 700 of file slepc_solver.h.

template<typename OutputVector >
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 738 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 116 of file slepc_solver.cc.

template<typename Vector >
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 794 of file slepc_solver.h.

void SLEPcWrappers::SolverBase::set_transformation ( SLEPcWrappers::TransformationBase this_transformation)

Set the spectral transformation to be used.

Definition at line 90 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 129 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.

Note
For other allowed values see the SLEPc documentation.

Definition at line 139 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.

Note
For other allowed values see the SLEPc documentation.

Definition at line 147 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 292 of file slepc_solver.cc.

SolverControl & SLEPcWrappers::SolverBase::control ( ) const

Access to the object that controls convergence.

Definition at line 318 of file slepc_solver.cc.

void SLEPcWrappers::SolverBase::solve ( const unsigned int  n_eigenpairs,
unsigned int *  n_converged 
)
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 154 of file slepc_solver.cc.

void SLEPcWrappers::SolverBase::get_eigenpair ( const unsigned int  index,
PetscScalar &  eigenvalues,
PETScWrappers::VectorBase eigenvectors 
)
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 248 of file slepc_solver.cc.

void SLEPcWrappers::SolverBase::get_eigenpair ( const unsigned int  index,
double &  real_eigenvalues,
double &  imag_eigenvalues,
PETScWrappers::VectorBase real_eigenvectors,
PETScWrappers::VectorBase imag_eigenvectors 
)
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 260 of file slepc_solver.cc.

void SLEPcWrappers::SolverBase::set_matrices ( const PETScWrappers::MatrixBase A)
protected

Initialize solver for the linear system \(Ax=\lambda x\). (Note: this is required before calling solve ())

Definition at line 73 of file slepc_solver.cc.

void SLEPcWrappers::SolverBase::set_matrices ( const PETScWrappers::MatrixBase A,
const PETScWrappers::MatrixBase B 
)
protected

Same as above, but here initialize solver for the linear system \(A x=\lambda B x\).

Definition at line 81 of file slepc_solver.cc.

int SLEPcWrappers::SolverBase::convergence_test ( EPS  eps,
PetscScalar  real_eigenvalue,
PetscScalar  imag_eigenvalue,
PetscReal  residual_norm,
PetscReal *  estimated_error,
void *  solver_control 
)
staticprivate

A function that can be used in SLEPc as a callback to check on convergence.

Note
This function is not used currently.

Definition at line 324 of file slepc_solver.cc.

Member Data Documentation

SolverControl& SLEPcWrappers::SolverBase::solver_control
protected

Reference to the object that controls convergence of the iterative solver.

Definition at line 304 of file slepc_solver.h.

const MPI_Comm SLEPcWrappers::SolverBase::mpi_communicator
protected

Copy of the MPI communicator object to be used for the solver.

Definition at line 309 of file slepc_solver.h.

EPS SLEPcWrappers::SolverBase::eps
protected

Objects for Eigenvalue Problem Solver.

Definition at line 361 of file slepc_solver.h.

EPSConvergedReason SLEPcWrappers::SolverBase::reason
private

Convergence reason.

Definition at line 367 of file slepc_solver.h.


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