Classes | |
struct | AdditionalData |
class | ExcTooFewTmpVectors |
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) |
Protected Member Functions | |
virtual double | criterion () |
void | givens_rotation (Vector< double > &h, Vector< double > &b, Vector< double > &ci, Vector< double > &si, int col) const |
Protected Attributes | |
AdditionalData | additional_data |
FullMatrix< double > | H |
FullMatrix< double > | H1 |
Private Member Functions | |
SolverGMRES (const SolverGMRES< VECTOR > &) |
You have to give the maximum number of temporary vectors to the constructor which are to be used to do the orthogonalization. If the number of iterations needed to solve the problem to the given criterion, an intermediate solution is computed and a restart is performed. If you don't want to use the restarted method, you can limit the number of iterations (stated in the SolverControl
object given to the constructor) to be below the number of temporary vectors minus three. Note the subtraction, which is due to the fact that three vectors are used for other purposes, so the number of iterations before a restart occurs is less by three than the total number of temporary vectors. If the size of the matrix is smaller than the maximum number of temporary vectors, then fewer vectors are allocated since GMRES can then be used as a direct solver.
Note that restarts don't compensate temporary vectors very well, i.e. giving too few temporary vectors will increase the necessary iteration steps greatly; it is not uncommon that you will need more iterations than there are degrees of freedom in your solution vector, if the number of temporary vectors is lower than the size of the vector, even though GMRES is an exact solver if a sufficient number of temporary vectors is given. You should therefore always give as many temporary vectors as you can, unless you are limited by the available memory; only then should you start to trade computational speed against memory. One of the few other possibilities is to use a good preconditioner.
Like all other solver classes, this class has a local structure called AdditionalData
which is used to pass additional parameters to the solver, like the number of temporary vectors for example. We use this additional structure instead of passing these values directly to the constructor because this makes the use of the SolverSelector
and other classes much easier and guarantees that these will continue to work even if number or type of the additional parameters for a certain solver changes.
For the GMRes method, the AdditionalData
structure contains the number of temporary vectors as commented upon above. By default, the number of these vectors is set to 30. The AdditionalData
also containes a flag indicating the use of 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. By default, the left preconditioned GMRes uses the preconditioned residual and the right preconditioned GMRes uses the normal, i.e. unpreconditioned, residual as stopping criterion. If the use_default_residual
flag is false
, the left preconditioned GMRes uses as stopping criterion the unpreconditioned residual and the right preconditioned GMRes the preconditioned residual. But be aware that the non-default residuals are not automatically computed by the GMRes method but need to be computed in addition. This (especially for the left preconditioned GMRes) might lead to a significant loss in the solver performance. Therefore, the user should set use_default_residual=false
only for debugging/testing purposes.
For the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class.
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.
SolverGMRES< VECTOR >::SolverGMRES | ( | const SolverGMRES< VECTOR > & | ) | [private] |
No copy constructor.
void SolverGMRES< VECTOR >::solve | ( | const MATRIX & | A, | |
VECTOR & | x, | |||
const VECTOR & | b, | |||
const PRECONDITIONER & | precondition | |||
) | [inline] |
Solve the linear system for x.
Referenced by SolverSelector< VECTOR >::solve(), and EigenInverse< VECTOR >::solve().
virtual double SolverGMRES< VECTOR >::criterion | ( | ) | [protected, virtual] |
Implementation of the computation of the norm of the residual.
void SolverGMRES< VECTOR >::givens_rotation | ( | Vector< double > & | h, | |
Vector< double > & | b, | |||
Vector< double > & | ci, | |||
Vector< double > & | si, | |||
int | col | |||
) | const [protected] |
Transformation of an upper Hessenberg matrix into tridiagonal structure by givens rotation of the last column
AdditionalData SolverGMRES< VECTOR >::additional_data [protected] |
Includes the maximum number of tmp vectors.
FullMatrix<double> SolverGMRES< VECTOR >::H [protected] |
Projected system matrix
FullMatrix<double> SolverGMRES< VECTOR >::H1 [protected] |
Auxiliary matrix for inverting H