SolverGMRES< VECTOR > Class Template Reference
[Linear solver classes]

Inheritance diagram for SolverGMRES< VECTOR >:

Inheritance graph
[legend]

List of all members.

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)

Static Public Member Functions

::ExceptionBaseExcTooFewTmpVectors (int arg1) throw (errortext << "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." )

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< doubleH
FullMatrix< doubleH1

Private Member Functions

 SolverGMRES (const SolverGMRES< VECTOR > &)

Classes

struct  AdditionalData


Detailed Description

template<class VECTOR = Vector<double>>
class 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. Addinitonally, it allows you to choose bet 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.

Left versus right preconditioning

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 size of the Arnoldi basis

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.

Author:
Wolfgang Bangerth, Guido Kanschat, Ralf Hartmann.

Constructor & Destructor Documentation

template<class VECTOR = Vector<double>>
SolverGMRES< VECTOR >::SolverGMRES ( SolverControl cn,
VectorMemory< VECTOR > &  mem,
const AdditionalData data = AdditionalData() 
)

Constructor.

template<class VECTOR = Vector<double>>
SolverGMRES< VECTOR >::SolverGMRES ( SolverControl cn,
const AdditionalData data = AdditionalData() 
)

Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.

template<class VECTOR = Vector<double>>
SolverGMRES< VECTOR >::SolverGMRES ( const SolverGMRES< VECTOR > &   )  [private]

No copy constructor.


Member Function Documentation

template<class VECTOR = Vector<double>>
template<class MATRIX, class PRECONDITIONER>
void SolverGMRES< VECTOR >::solve ( const MATRIX &  A,
VECTOR &  x,
const VECTOR &  b,
const PRECONDITIONER &  precondition 
) [inline]

Solve the linear system $Ax=b$ for x.

Referenced by SolverSelector< VECTOR >::solve(), and EigenInverse< VECTOR >::solve().

template<class VECTOR = Vector<double>>
::ExceptionBase& SolverGMRES< VECTOR >::ExcTooFewTmpVectors ( int  arg1  )  throw (errortext << "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." ) [static]

template<class VECTOR = Vector<double>>
virtual double SolverGMRES< VECTOR >::criterion (  )  [protected, virtual]

Implementation of the computation of the norm of the residual.

template<class VECTOR = Vector<double>>
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


Member Data Documentation

template<class VECTOR = Vector<double>>
AdditionalData SolverGMRES< VECTOR >::additional_data [protected]

Includes the maximum number of tmp vectors.

template<class VECTOR = Vector<double>>
FullMatrix<double> SolverGMRES< VECTOR >::H [protected]

Projected system matrix

template<class VECTOR = Vector<double>>
FullMatrix<double> SolverGMRES< VECTOR >::H1 [protected]

Auxiliary matrix for inverting H


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

deal.II documentation generated on Sun Sep 5 23:06:43 2010 by doxygen 1.5.6