Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/solver_gmres.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
SolverFGMRES (SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData()) | |
SolverFGMRES (SolverControl &cn, const AdditionalData &data=AdditionalData()) | |
template<typename MatrixType , typename PreconditionerType > | |
void | solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner) |
Public Member Functions inherited from Solver< VectorType > | |
Solver (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory) | |
Solver (SolverControl &solver_control) | |
boost::signals2::connection | connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate)> &slot) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
AdditionalData | additional_data |
FullMatrix< double > | H |
FullMatrix< double > | H1 |
Additional Inherited Members | |
Public Types inherited from Solver< VectorType > | |
typedef VectorType | vector_type |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Attributes inherited from Solver< VectorType > | |
GrowingVectorMemory< VectorType > | static_vector_memory |
VectorMemory< VectorType > & | memory |
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > | iteration_status |
Implementation of the Generalized minimal residual method with flexible preconditioning method.
This version of the GMRES method allows for the use of a different preconditioner in each iteration step. Therefore, it is also more robust with respect to inaccurate evaluation of the preconditioner. An important application is also the use of a Krylov space method inside the preconditioner. As opposed to SolverGMRES which allows one to choose between left and right preconditioning, this solver always applies the preconditioner from the right.
FGMRES needs two vectors in each iteration steps yielding a total of 2*SolverFGMRES::AdditionalData::max_basis_size+1
auxiliary vectors.
Caveat: Documentation of this class is not up to date. There are also a few parameters of GMRES we would like to introduce here.
Definition at line 452 of file solver_gmres.h.
SolverFGMRES< VectorType >::SolverFGMRES | ( | SolverControl & | cn, |
VectorMemory< VectorType > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverFGMRES< VectorType >::SolverFGMRES | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
void SolverFGMRES< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve the linear system \(Ax=b\) for x.
|
private |
Additional flags.
Definition at line 505 of file solver_gmres.h.
|
private |
Projected system matrix
Definition at line 510 of file solver_gmres.h.
|
private |
Auxiliary matrix for inverting H
Definition at line 515 of file solver_gmres.h.