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 Index &iters,
const Index &restart,
typename Dest::RealScalar & tol_error) {
62 typedef typename Dest::RealScalar RealScalar;
63 typedef typename Dest::Scalar Scalar;
64 typedef Matrix < Scalar, Dynamic, 1 > VectorType;
65 typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
67 RealScalar tol = tol_error;
68 const Index maxIters = iters;
71 const Index m = mat.rows();
74 VectorType p0 = rhs - mat*x;
75 VectorType r0 = precond.solve(p0);
77 const RealScalar r0Norm = r0.norm();
87 FMatrixType H = FMatrixType::Zero(m, restart + 1);
88 VectorType w = VectorType::Zero(restart + 1);
89 VectorType tau = VectorType::Zero(restart + 1);
92 std::vector < JacobiRotation < Scalar > > G(restart);
95 VectorType t(m), v(m), workspace(m), x_new(m);
98 Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
100 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
103 for (Index k = 1; k <= restart; ++k)
107 v = VectorType::Unit(m, k - 1);
111 for (Index i = k - 1; i >= 0; --i) {
112 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
116 t.noalias() = mat * v;
117 v = precond.solve(t);
121 for (Index i = 0; i < k; ++i) {
122 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
125 if (v.tail(m - k).norm() != 0.0)
130 Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
131 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
134 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
140 for (Index i = 0; i < k - 1; ++i)
143 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
147 if (k<m && v(k) != (Scalar) 0)
150 G[k - 1].makeGivens(v(k - 1), v(k));
153 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
154 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
158 H.col(k-1).head(k) = v.head(k);
160 bool stop = (k==m || abs(w(k)) < tol * r0Norm || iters == maxIters);
162 if (stop || k == restart)
165 Ref<VectorType> y = w.head(k);
166 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
170 for (Index i = k - 1; i >= 0; --i)
174 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
188 p0.noalias() = rhs - mat*x;
189 r0 = precond.solve(p0);
197 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
209 template<
typename _MatrixType,
210 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
215 template<
typename _MatrixType,
typename _Preconditioner>
216 struct traits<
GMRES<_MatrixType,_Preconditioner> >
218 typedef _MatrixType MatrixType;
219 typedef _Preconditioner Preconditioner;
258 template<
typename _MatrixType,
typename _Preconditioner>
259 class GMRES :
public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
261 typedef IterativeSolverBase<GMRES> Base;
264 using Base::m_iterations;
266 using Base::m_isInitialized;
272 using Base::_solve_impl;
273 typedef _MatrixType MatrixType;
274 typedef typename MatrixType::Scalar Scalar;
275 typedef typename MatrixType::RealScalar RealScalar;
276 typedef _Preconditioner Preconditioner;
293 template<
typename MatrixDerived>
294 explicit GMRES(
const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
308 template<
typename Rhs,
typename Dest>
309 void _solve_with_guess_impl(
const Rhs& b, Dest& x)
const
312 for(Index j=0; j<b.cols(); ++j)
314 m_iterations = Base::maxIterations();
315 m_error = Base::m_tolerance;
317 typename Dest::ColXpr xj(x,j);
318 if(!internal::gmres(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
321 m_info = failed ? NumericalIssue
322 : m_error <= Base::m_tolerance ? Success
324 m_isInitialized =
true;
328 template<
typename Rhs,
typename Dest>
329 void _solve_impl(
const Rhs& b, MatrixBase<Dest> &x)
const
332 if(x.squaredNorm() == 0)
return;
333 _solve_with_guess_impl(b,x.derived());
342 #endif // EIGEN_GMRES_H
GMRES()
Definition: GMRES.h:281
Index get_restart()
Definition: GMRES.h:300
Namespace containing all symbols from the Eigen library.
Definition: CXX11Meta.h:13
GMRES(const EigenBase< MatrixDerived > &A)
Definition: GMRES.h:294
void set_restart(const Index restart)
Definition: GMRES.h:305
A GMRES solver for sparse square problems.
Definition: GMRES.h:211