Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Types
Eigen::DGMRES< _MatrixType, _Preconditioner > Class Template Reference

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. More...

#include <DGMRES.h>

Inheritance diagram for Eigen::DGMRES< _MatrixType, _Preconditioner >:
Inheritance graph
[legend]

List of all members.

Public Types

typedef Matrix< std::complex
< RealScalar >, Dynamic, 1 > 
ComplexVector
typedef Matrix< Scalar,
Dynamic, Dynamic
DenseMatrix
typedef Matrix< RealScalar,
Dynamic, Dynamic
DenseRealMatrix
typedef Matrix< RealScalar,
Dynamic, 1 > 
DenseRealVector
typedef Matrix< Scalar,
Dynamic, 1 > 
DenseVector
typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef _Preconditioner Preconditioner
typedef MatrixType::RealScalar RealScalar
typedef MatrixType::Scalar Scalar

Public Member Functions

template<typename Rhs , typename Dest >
void _solve (const Rhs &b, Dest &x) const
template<typename Rhs , typename Dest >
void _solveWithGuess (const Rhs &b, Dest &x) const
int deflSize ()
 DGMRES ()
 DGMRES (const MatrixType &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
 ~DGMRES ()

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,.
template<typename RhsType , typename DestType >
int dgmresApplyDeflation (const RhsType &In, DestType &Out) const
int dgmresComputeDeflationData (const MatrixType &mat, const Preconditioner &precond, const Index &it, Index &neig) const
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.
void dgmresInitDeflation (Index &rows) const
ComplexVector schurValues (const ComplexSchur< DenseMatrix > &schurofH) const
ComplexVector schurValues (const RealSchur< DenseMatrix > &schurofH) const

Protected Attributes

bool m_force
DenseMatrix m_H
DenseMatrix m_Hes
bool m_isDeflAllocated
bool m_isDeflInitialized
RealScalar m_lambdaN
PartialPivLU< DenseMatrixm_luT
int m_maxNeig
DenseMatrix m_MU
int m_neig
int m_r
Index m_restart
RealScalar m_smv
DenseMatrix m_T
DenseMatrix m_U
DenseMatrix m_V

Private Types

typedef IterativeSolverBase
< DGMRES
Base

Detailed Description

template<typename _MatrixType, typename _Preconditioner>
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.

Definition at line 101 of file DGMRES.h.


Member Typedef Documentation

template<typename _MatrixType , typename _Preconditioner >
typedef IterativeSolverBase<DGMRES> Eigen::DGMRES< _MatrixType, _Preconditioner >::Base [private]

Definition at line 103 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef Matrix<std::complex<RealScalar>, Dynamic, 1> Eigen::DGMRES< _MatrixType, _Preconditioner >::ComplexVector

Definition at line 120 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef Matrix<Scalar,Dynamic,Dynamic> Eigen::DGMRES< _MatrixType, _Preconditioner >::DenseMatrix

Definition at line 116 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef Matrix<RealScalar,Dynamic,Dynamic> Eigen::DGMRES< _MatrixType, _Preconditioner >::DenseRealMatrix

Definition at line 117 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef Matrix<RealScalar,Dynamic,1> Eigen::DGMRES< _MatrixType, _Preconditioner >::DenseRealVector

Definition at line 119 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef Matrix<Scalar,Dynamic,1> Eigen::DGMRES< _MatrixType, _Preconditioner >::DenseVector

Definition at line 118 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::Index Eigen::DGMRES< _MatrixType, _Preconditioner >::Index
template<typename _MatrixType , typename _Preconditioner >
typedef _MatrixType Eigen::DGMRES< _MatrixType, _Preconditioner >::MatrixType
template<typename _MatrixType , typename _Preconditioner >
typedef _Preconditioner Eigen::DGMRES< _MatrixType, _Preconditioner >::Preconditioner
template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::RealScalar Eigen::DGMRES< _MatrixType, _Preconditioner >::RealScalar
template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::Scalar Eigen::DGMRES< _MatrixType, _Preconditioner >::Scalar

Constructor & Destructor Documentation

template<typename _MatrixType , typename _Preconditioner >
Eigen::DGMRES< _MatrixType, _Preconditioner >::DGMRES ( ) [inline]

Default constructor.

Definition at line 124 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
Eigen::DGMRES< _MatrixType, _Preconditioner >::DGMRES ( const MatrixType A) [inline]

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.

Definition at line 136 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
Eigen::DGMRES< _MatrixType, _Preconditioner >::~DGMRES ( ) [inline]

Definition at line 139 of file DGMRES.h.


Member Function Documentation

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::_solve ( const Rhs &  b,
Dest &  x 
) const [inline]

Definition at line 178 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::_solveWithGuess ( const Rhs &  b,
Dest &  x 
) const [inline]

Definition at line 159 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::deflSize ( ) [inline]

Get the size of the deflation subspace size

Definition at line 205 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::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.

Definition at line 256 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename RhsType , typename DestType >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::dgmresApplyDeflation ( const RhsType &  In,
DestType &  Out 
) const [protected]

Definition at line 518 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::dgmresComputeDeflationData ( const MatrixType mat,
const Preconditioner precond,
const Index it,
Index neig 
) const [protected]

Definition at line 434 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Dest >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::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

Definition at line 300 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::dgmresInitDeflation ( Index rows) const [protected]

Definition at line 392 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::restart ( ) [inline]

Get the restart value

Definition at line 186 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DGMRES< _MatrixType, _Preconditioner >::ComplexVector Eigen::DGMRES< _MatrixType, _Preconditioner >::schurValues ( const ComplexSchur< DenseMatrix > &  schurofH) const [inline, protected]

Definition at line 402 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DGMRES< _MatrixType, _Preconditioner >::ComplexVector Eigen::DGMRES< _MatrixType, _Preconditioner >::schurValues ( const RealSchur< DenseMatrix > &  schurofH) const [inline, protected]

Definition at line 408 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::set_restart ( const int  restart) [inline]

Set the restart value (default is 30)

Definition at line 191 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::setEigenv ( const int  neig) [inline]

Set the number of eigenvalues to deflate at each restart

Definition at line 196 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
void Eigen::DGMRES< _MatrixType, _Preconditioner >::setMaxEigenv ( const int  maxNeig) [inline]

Set the maximum size of the deflation subspace

Definition at line 210 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Guess >
const internal::solve_retval_with_guess<DGMRES, Rhs, Guess> Eigen::DGMRES< _MatrixType, _Preconditioner >::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()

Definition at line 148 of file DGMRES.h.


Member Data Documentation

template<typename _MatrixType , typename _Preconditioner >
bool Eigen::DGMRES< _MatrixType, _Preconditioner >::m_force [mutable, protected]

Definition at line 245 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_H [mutable, protected]

Definition at line 229 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_Hes [mutable, protected]

Definition at line 230 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
bool Eigen::DGMRES< _MatrixType, _Preconditioner >::m_isDeflAllocated [mutable, protected]

Definition at line 240 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
bool Eigen::DGMRES< _MatrixType, _Preconditioner >::m_isDeflInitialized [mutable, protected]

Definition at line 241 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
RealScalar Eigen::DGMRES< _MatrixType, _Preconditioner >::m_lambdaN [mutable, protected]

Definition at line 239 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
PartialPivLU<DenseMatrix> Eigen::DGMRES< _MatrixType, _Preconditioner >::m_luT [mutable, protected]

Definition at line 235 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::m_maxNeig [mutable, protected]

Definition at line 238 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_MU [mutable, protected]

Definition at line 233 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::m_neig [mutable, protected]

Definition at line 236 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::DGMRES< _MatrixType, _Preconditioner >::m_r [mutable, protected]

Definition at line 237 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
Index Eigen::DGMRES< _MatrixType, _Preconditioner >::m_restart [mutable, protected]

Definition at line 231 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
RealScalar Eigen::DGMRES< _MatrixType, _Preconditioner >::m_smv [mutable, protected]

Definition at line 244 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_T [mutable, protected]

Definition at line 234 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_U [mutable, protected]

Definition at line 232 of file DGMRES.h.

template<typename _MatrixType , typename _Preconditioner >
DenseMatrix Eigen::DGMRES< _MatrixType, _Preconditioner >::m_V [mutable, protected]

Definition at line 228 of file DGMRES.h.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 12:01:51