![]() |
Reference documentation for deal.II version 8.1.0
|
#include <solver_gmres.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
SolverGMRES (SolverControl &cn, VectorMemory< VECTOR > &mem, const AdditionalData &data=AdditionalData()) | |
SolverGMRES (SolverControl &cn, const AdditionalData &data=AdditionalData()) | |
template<class MATRIX , class PRECONDITIONER > | |
void | solve (const MATRIX &A, VECTOR &x, const VECTOR &b, const PRECONDITIONER &precondition) |
DeclException1 (ExcTooFewTmpVectors, int,<< "The number of temporary vectors you gave ("<< arg1<< ") is too small. It should be at least 10 for "<< "any results, and much more for reasonable ones.") | |
![]() | |
Solver (SolverControl &solver_control, VectorMemory< VECTOR > &vector_memory) | |
Solver (SolverControl &solver_control) | |
SolverControl & | control () const |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
virtual double | criterion () |
void | givens_rotation (Vector< double > &h, Vector< double > &b, Vector< double > &ci, Vector< double > &si, int col) const |
Static Protected Member Functions | |
static double | modified_gram_schmidt (const internal::SolverGMRES::TmpVectors< VECTOR > &orthogonal_vectors, const unsigned int dim, const unsigned int accumulated_iterations, VECTOR &vv, Vector< double > &h, bool &re_orthogonalize) |
Protected Attributes | |
AdditionalData | additional_data |
FullMatrix< double > | H |
FullMatrix< double > | H1 |
![]() | |
GrowingVectorMemory< VECTOR > | static_vector_memory |
SolverControl & | cntrl |
VectorMemory< VECTOR > & | memory |
Private Member Functions | |
SolverGMRES (const SolverGMRES< VECTOR > &) | |
Implementation of the Restarted Preconditioned Direct Generalized Minimal Residual Method. The stopping criterion is the norm of the residual.
The AdditionalData structure contains the number of temporary vectors used. The size of the Arnoldi basis is this number minus three. Additionally, it allows you to choose between right or left preconditioning. The default is left preconditioning. Finally it includes a flag indicating whether or not the default residual is used as stopping criterion.
AdditionalData
allows you to choose between left and right preconditioning. As expected, this switches between solving for the systems P-1A and AP-1, respectively.
A second consequence is the type of residual which is used to measure convergence. With left preconditioning, this is the preconditioned residual, while with right preconditioning, it is the residual of the unpreconditioned system.
Optionally, this behavior can be overridden by using the flag AdditionalData::use_default_residual. A true
value refers to the behavior described in the previous paragraph, while false
reverts it. Be aware though that additional residuals have to be computed in this case, impeding the overall performance of the solver.
The maximal basis size is controlled by AdditionalData::max_n_tmp_vectors, and it is this number minus 2. If the number of iteration steps exceeds this number, all basis vectors are discarded and the iteration starts anew from the approximation obtained so far.
Note that the minimizing property of GMRes only pertains to the Krylov space spanned by the Arnoldi basis. Therefore, restarted GMRes is not minimizing anymore. The choice of the basis length is a trade-off between memory consumption and convergence speed, since a longer basis means minimization over a larger space.
For the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class.
Definition at line 156 of file solver_gmres.h.
SolverGMRES< VECTOR >::SolverGMRES | ( | SolverControl & | cn, |
VectorMemory< VECTOR > & | mem, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverGMRES< VECTOR >::SolverGMRES | ( | SolverControl & | cn, |
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
|
private |
No copy constructor.
void SolverGMRES< VECTOR >::solve | ( | const MATRIX & | A, |
VECTOR & | x, | ||
const VECTOR & | b, | ||
const PRECONDITIONER & | precondition | ||
) |
Solve the linear system for x.
|
protectedvirtual |
Implementation of the computation of the norm of the residual.
|
protected |
Transformation of an upper Hessenberg matrix into tridiagonal structure by givens rotation of the last column
|
staticprotected |
Orthogonalize the vector vv
against the dim
(orthogonal) vectors given by the first argument using the modified Gram-Schmidt algorithm. The factors used for orthogonalization are stored in h
. The boolean re_orthogonalize
specifies whether the modified Gram-Schmidt algorithm should be applied twice. The algorithm checks loss of orthogonality in the procedure every fifth step and sets the flag to true in that case. All subsequent iterations use re-orthogonalization.
|
protected |
Includes the maximum number of tmp vectors.
Definition at line 239 of file solver_gmres.h.
|
protected |
Projected system matrix
Definition at line 274 of file solver_gmres.h.
|
protected |
Auxiliary matrix for inverting H
Definition at line 279 of file solver_gmres.h.