Reference documentation for deal.II version GIT relicensing-216-gcda843ca70 2024-03-28 12:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | Private Attributes | List of all members
SLEPcWrappers::SolverArnoldi Class Reference

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

Inheritance diagram for SLEPcWrappers::SolverArnoldi:
Inheritance graph
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 SolverArnoldi (SolverControl &cn, const MPI_Comm mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
 
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)
 
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 ::ExceptionBaseExcSLEPcWrappersUsageError ()
 
static ::ExceptionBaseExcSLEPcError (int arg1)
 
static ::ExceptionBaseExcSLEPcEigenvectorConvergenceMismatchError (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

const AdditionalData additional_data
 
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

An implementation of the solver interface using the SLEPc Arnoldi solver. Usage: All spectrum, all problem types, complex.

For examples of how this and its sibling classes can be used, including how to provide preconditioners to the matrix of which eigenvalues are to be computed, see the documentation of the SolverBase class as well as the extensive discussions in the documentation of the SLEPcWrappers namespace.

Definition at line 432 of file slepc_solver.h.

Constructor & Destructor Documentation

◆ SolverArnoldi()

SLEPcWrappers::SolverArnoldi::SolverArnoldi ( SolverControl cn,
const MPI_Comm  mpi_communicator = PETSC_COMM_SELF,
const AdditionalData data = AdditionalData() 
)
explicit

SLEPc solvers will want to have an MPI communicator context over which computations are parallelized. By default, this carries the same behavior as the PETScWrappers, but you can change that.

Definition at line 382 of file slepc_solver.cc.

Member Function Documentation

◆ solve() [1/4]

template<typename OutputVector >
void SolverBase< OutputVector >::solve ( const PETScWrappers::MatrixBase A,
std::vector< PetscScalar > &  eigenvalues,
std::vector< OutputVector > &  eigenvectors,
const unsigned int  n_eigenpairs = 1 
)
inherited

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

◆ solve() [2/4]

template<typename OutputVector >
void SolverBase< OutputVector >::solve ( const PETScWrappers::MatrixBase A,
const PETScWrappers::MatrixBase B,
std::vector< PetscScalar > &  eigenvalues,
std::vector< OutputVector > &  eigenvectors,
const unsigned int  n_eigenpairs = 1 
)
inherited

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

◆ solve() [3/4]

template<typename OutputVector >
void SolverBase< OutputVector >::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 
)
inherited

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

◆ solve() [4/4]

void SolverBase< VectorType >::solve ( const unsigned int  n_eigenpairs,
unsigned int n_converged 
)
protectedinherited

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

◆ set_initial_space()

template<typename Vector >
void SolverBase< Vector >::set_initial_space ( const std::vector< Vector > &  initial_space)
inherited

Set the initial vector space for the solver.

By default, SLEPc initializes the starting vector or the initial subspace randomly.

Definition at line 831 of file slepc_solver.h.

◆ set_transformation()

void SolverBase< VectorType >::set_transformation ( SLEPcWrappers::TransformationBase this_transformation)
inherited

Set the spectral transformation to be used.

Definition at line 97 of file slepc_solver.cc.

◆ set_target_eigenvalue()

void SolverBase< VectorType >::set_target_eigenvalue ( const PetscScalar &  this_target)
inherited

Set target eigenvalues in the spectrum to be computed. By default, no target is set.

Definition at line 125 of file slepc_solver.cc.

◆ set_which_eigenpairs()

void SolverBase< VectorType >::set_which_eigenpairs ( EPSWhich  set_which)
inherited

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

◆ set_problem_type()

void SolverBase< VectorType >::set_problem_type ( EPSProblemType  set_problem)
inherited

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.

◆ get_solver_state()

void SolverBase< VectorType >::get_solver_state ( const SolverControl::State  state)
inherited

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

◆ control()

SolverControl & SolverBase< VectorType >::control ( ) const
inherited

Access to the object that controls convergence.

Definition at line 334 of file slepc_solver.cc.

◆ get_eigenpair() [1/2]

void SolverBase< VectorType >::get_eigenpair ( const unsigned int  index,
PetscScalar &  eigenvalues,
PETScWrappers::VectorBase eigenvectors 
)
protectedinherited

Access the real parts of solutions for a solved eigenvector problem, pair index solutions, \(\text{index}\,\in\,0\dots \mathrm{n\_converged}-1\).

Definition at line 259 of file slepc_solver.cc.

◆ get_eigenpair() [2/2]

void SolverBase< VectorType >::get_eigenpair ( const unsigned int  index,
double &  real_eigenvalues,
double &  imag_eigenvalues,
PETScWrappers::VectorBase real_eigenvectors,
PETScWrappers::VectorBase imag_eigenvectors 
)
protectedinherited

Access the real and imaginary parts of solutions for a solved eigenvector problem, pair index solutions, \(\text{index}\,\in\,0\dots \mathrm{n\_converged}-1\).

Definition at line 272 of file slepc_solver.cc.

◆ set_matrices() [1/2]

void SolverBase< VectorType >::set_matrices ( const PETScWrappers::MatrixBase A)
protectedinherited

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

Definition at line 76 of file slepc_solver.cc.

◆ set_matrices() [2/2]

void SolverBase< VectorType >::set_matrices ( const PETScWrappers::MatrixBase A,
const PETScWrappers::MatrixBase B 
)
protectedinherited

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

Definition at line 86 of file slepc_solver.cc.

◆ convergence_test()

int SolverBase< VectorType >::convergence_test ( EPS  eps,
PetscScalar  real_eigenvalue,
PetscScalar  imag_eigenvalue,
PetscReal  residual_norm,
PetscReal *  estimated_error,
void *  solver_control 
)
staticprivateinherited

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

Member Data Documentation

◆ additional_data

const AdditionalData SLEPcWrappers::SolverArnoldi::additional_data
protected

Store a copy of the flags for this particular solver.

Definition at line 466 of file slepc_solver.h.

◆ solver_control

SolverControl& SLEPcWrappers::SolverBase::solver_control
protectedinherited

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

Definition at line 296 of file slepc_solver.h.

◆ mpi_communicator

const MPI_Comm SLEPcWrappers::SolverBase::mpi_communicator
protectedinherited

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

Definition at line 301 of file slepc_solver.h.

◆ eps

EPS SLEPcWrappers::SolverBase::eps
protectedinherited

Objects for Eigenvalue Problem Solver.

Definition at line 353 of file slepc_solver.h.

◆ reason

EPSConvergedReason SLEPcWrappers::SolverBase::reason
privateinherited

Convergence reason.

Definition at line 359 of file slepc_solver.h.


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