55 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
56 bool gmres(
const MatrixType & mat,
const Rhs & rhs, Dest & x,
const Preconditioner & precond,
57 int &iters,
const int &restart,
typename Dest::RealScalar & tol_error) {
62 typedef typename Dest::RealScalar RealScalar;
67 RealScalar tol = tol_error;
68 const int maxIters = iters;
71 const int m = mat.rows();
73 VectorType p0 = rhs - mat*x;
74 VectorType r0 = precond.solve(p0);
77 if(abs(r0.norm()) < tol) {
81 VectorType w = VectorType::Zero(restart + 1);
83 FMatrixType H = FMatrixType::Zero(m, restart + 1);
84 VectorType tau = VectorType::Zero(restart + 1);
85 std::vector < JacobiRotation < Scalar > > G(restart);
90 r0.makeHouseholder(e, tau.coeffRef(0), beta);
92 H.bottomLeftCorner(m - 1, 1) = e;
94 for (
int k = 1; k <= restart; ++k) {
98 VectorType v = VectorType::Unit(m, k - 1), workspace(m);
101 for (
int i = k - 1; i >= 0; --i) {
102 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
110 for (
int i = 0; i < k; ++i) {
111 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
114 if (v.tail(m - k).norm() != 0.0) {
119 VectorType e(m - k - 1);
121 v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
122 H.col(k).tail(m - k - 1) = e;
125 v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
131 for (
int i = 0; i < k - 1; ++i) {
133 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
137 if (k<m && v(k) != (Scalar) 0) {
139 G[k - 1].makeGivens(v(k - 1), v(k));
142 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
143 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
148 H.col(k - 1).head(k) = v.head(k);
150 bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
152 if (stop || k == restart) {
155 VectorType y = w.head(k);
156 H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
159 VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
162 x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
164 for (
int i = k - 2; i >= 0; --i) {
165 x_new += y(i) * VectorType::Unit(m, i);
167 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
179 VectorType p1=precond.solve(p0);
182 w = VectorType::Zero(restart + 1);
183 H = FMatrixType::Zero(m, restart + 1);
184 tau = VectorType::Zero(restart + 1);
188 r0.makeHouseholder(e, tau.coeffRef(0), beta);
190 H.bottomLeftCorner(m - 1, 1) = e;
206 template<
typename _MatrixType,
212 template<
typename _MatrixType,
typename _Preconditioner>
215 typedef _MatrixType MatrixType;
216 typedef _Preconditioner Preconditioner;
253 template<
typename _MatrixType,
typename _Preconditioner>
257 using Base::mp_matrix;
259 using Base::m_iterations;
261 using Base::m_isInitialized;
267 typedef _MatrixType MatrixType;
269 typedef typename MatrixType::Index Index;
270 typedef typename MatrixType::RealScalar RealScalar;
271 typedef _Preconditioner Preconditioner;
288 template<
typename MatrixDerived>
307 template<
typename Rhs,
typename Guess>
311 eigen_assert(m_isInitialized &&
"GMRES is not initialized.");
312 eigen_assert(Base::rows()==b.rows()
313 &&
"GMRES::solve(): invalid number of rows of the right hand side matrix b");
315 <
GMRES, Rhs, Guess>(*
this, b.derived(), x0);
319 template<
typename Rhs,
typename Dest>
320 void _solveWithGuess(
const Rhs& b, Dest& x)
const 323 for(
int j=0; j<b.cols(); ++j)
325 m_iterations = Base::maxIterations();
326 m_error = Base::m_tolerance;
328 typename Dest::ColXpr xj(x,j);
329 if(!
internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
333 : m_error <= Base::m_tolerance ?
Success 335 m_isInitialized =
true;
339 template<
typename Rhs,
typename Dest>
340 void _solve(
const Rhs& b, Dest& x)
const 343 if(x.squaredNorm() == 0)
return;
344 _solveWithGuess(b,x);
354 template<
typename _MatrixType,
typename _Preconditioner,
typename Rhs>
359 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
361 template<
typename Dest>
void evalTo(Dest& dst)
const 363 dec()._solve(rhs(),dst);
371 #endif // EIGEN_GMRES_H A preconditioner based on the digonal entries.
Definition: BasicPreconditioners.h:33
Definition: gtest_unittest.cc:5031
Definition: ForwardDeclarations.h:124
GMRES()
Default constructor.
Definition: GMRES.h:276
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: TestIMU_Common.h:87
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, int &iters, const int &restart, typename Dest::RealScalar &tol_error)
Generalized Minimal Residual Algorithm based on the Arnoldi algorithm implemented with Householder re...
Definition: GMRES.h:56
const internal::solve_retval_with_guess< GMRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: GMRES.h:309
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:26
int get_restart()
Get the number of iterations after that a restart is performed.
Definition: GMRES.h:295
The provided data did not satisfy the prerequisites.
Definition: Constants.h:378
Definition: SparseSolve.h:86
GMRES(const EigenBase< MatrixDerived > &A)
Initialize the solver with matrix A for further Ax=b solving.
Definition: GMRES.h:289
Computation was successful.
Definition: Constants.h:376
Definition: BandTriangularSolver.h:13
A GMRES solver for sparse square problems.
Definition: GMRES.h:208
Definition: ForwardDeclarations.h:125
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:21
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Iterative procedure did not converge.
Definition: Constants.h:380
Definition: ForwardDeclarations.h:17
void set_restart(const int restart)
Set the number of iterations after that a restart is performed.
Definition: GMRES.h:300
double Scalar
Common scalar type.
Definition: FlexibleKalmanBase.h:48