All Classes Namespaces Functions Variables Typedefs Enumerator Groups Pages
DGMRES< _MatrixType, _Preconditioner > Class Template Reference

Detailed Description

template<typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
class Eigen::DGMRES< _MatrixType, _Preconditioner >

A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse linear systems. The basis is built with modified Gram-Schmidt. At each restart, a few approximated eigenvectors corresponding to the smallest eigenvalues are used to build a preconditioner for the next cycle. This preconditioner for deflation can be combined with any other preconditioner, the IncompleteLUT for instance. The preconditioner is applied at right of the matrix and the combination is multiplicative.

Template Parameters
_MatrixTypethe type of the sparse matrix A, can be a dense or a sparse matrix.
_Preconditionerthe type of the preconditioner. Default is DiagonalPreconditioner Typical usage :
* SparseMatrix<double> A;
* VectorXd x, b;
* //Fill A and b ...
* DGMRES<SparseMatrix<double> > solver;
* solver.set_restart(30); // Set restarting value
* solver.setEigenv(1); // Set the number of eigenvalues to deflate
* solver.compute(A);
* x = solver.solve(b);
*

References : [1] D. NUENTSA WAKAM and F. PACULL, Memory Efficient Hybrid Algebraic Solvers for Linear Systems Arising from Compressible Flows, Computers and Fluids, In Press, http://dx.doi.org/10.1016/j.compfluid.2012.03.023 [2] K. Burrage and J. Erhel, On the performance of various adaptive preconditioned GMRES strategies, 5(1998), 101-121. [3] J. Erhel, K. Burrage and B. Pohl, Restarted GMRES preconditioned by deflation,J. Computational and Applied Mathematics, 69(1996), 303-318.

Inherits IterativeSolverBase< DGMRES< _MatrixType, _Preconditioner > >.

Public Member Functions

int deflSize ()
 
 DGMRES ()
 
template<typename MatrixDerived >
 DGMRES (const EigenBase< MatrixDerived > &A)
 
int restart ()
 
void set_restart (const int restart)
 
void setEigenv (const int neig)
 
void setMaxEigenv (const int maxNeig)
 
template<typename Rhs , typename Guess >
const
internal::solve_retval_with_guess
< DGMRES, Rhs, Guess > 
solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 

Protected Member Functions

template<typename Rhs , typename Dest >
void dgmres (const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond) const
 Perform several cycles of restarted GMRES with modified Gram Schmidt,. More...
 
template<typename Dest >
int dgmresCycle (const MatrixType &mat, const Preconditioner &precond, Dest &x, DenseVector &r0, RealScalar &beta, const RealScalar &normRhs, int &nbIts) const
 Perform one restart cycle of DGMRES. More...
 

Constructor & Destructor Documentation

DGMRES ( )
inline

Default constructor.

Referenced by DGMRES< _MatrixType, _Preconditioner >::solveWithGuess().

DGMRES ( const EigenBase< MatrixDerived > &  A)
inlineexplicit

Initialize the solver with matrix A for further Ax=b solving.

This constructor is a shortcut for the default constructor followed by a call to compute().

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

Member Function Documentation

int deflSize ( )
inline

Get the size of the deflation subspace size

void dgmres ( const MatrixType &  mat,
const Rhs &  rhs,
Dest &  x,
const Preconditioner &  precond 
) const
protected

Perform several cycles of restarted GMRES with modified Gram Schmidt,.

A right preconditioner is used combined with deflation.

int dgmresCycle ( const MatrixType &  mat,
const Preconditioner &  precond,
Dest &  x,
DenseVector &  r0,
RealScalar &  beta,
const RealScalar &  normRhs,
int &  nbIts 
) const
protected

Perform one restart cycle of DGMRES.

Parameters
matThe coefficient matrix
precondThe preconditioner
xthe new approximated solution
r0The initial residual vector
betaThe norm of the residual computed so far
normRhsThe norm of the right hand side vector
nbItsThe number of iterations
int restart ( )
inline

Get the restart value

Referenced by DGMRES< _MatrixType, _Preconditioner >::set_restart().

void set_restart ( const int  restart)
inline

Set the restart value (default is 30)

References DGMRES< _MatrixType, _Preconditioner >::restart().

void setEigenv ( const int  neig)
inline

Set the number of eigenvalues to deflate at each restart

void setMaxEigenv ( const int  maxNeig)
inline

Set the maximum size of the deflation subspace

const internal::solve_retval_with_guess<DGMRES, Rhs, Guess> solveWithGuess ( const MatrixBase< Rhs > &  b,
const Guess &  x0 
) const
inline
Returns
the solution x of $ A x = b $ using the current decomposition of A x0 as an initial solution.
See Also
compute()

References DGMRES< _MatrixType, _Preconditioner >::DGMRES().


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