14 #include <Eigen/IterativeLinearSolvers>
59 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
60 bool gmres(
const MatrixType & mat,
const Rhs & rhs, Dest & x,
const Preconditioner & precond,
61 Index &iters,
const Index &restart,
typename Dest::RealScalar & tol_error) {
66 typedef typename Dest::RealScalar RealScalar;
67 typedef typename Dest::Scalar Scalar;
68 typedef Matrix < Scalar, Dynamic, 1 > VectorType;
69 typedef Matrix < Scalar, Dynamic, Dynamic, ColMajor> FMatrixType;
71 RealScalar tol = tol_error;
72 const Index maxIters = iters;
75 const Index m = mat.rows();
78 VectorType p0 = rhs - mat*x;
79 VectorType r0 = precond.solve(p0);
81 const RealScalar r0Norm = r0.norm();
91 FMatrixType
H = FMatrixType::Zero(m, restart + 1);
92 VectorType w = VectorType::Zero(restart + 1);
93 VectorType tau = VectorType::Zero(restart + 1);
96 std::vector < JacobiRotation < Scalar > > G(restart);
99 VectorType t(m), v(m), workspace(m), x_new(m);
102 Ref<VectorType> H0_tail =
H.col(0).tail(m - 1);
104 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
107 for (Index k = 1; k <= restart; ++k)
111 v = VectorType::Unit(m, k - 1);
115 for (Index i = k - 1; i >= 0; --i) {
116 v.tail(m - i).applyHouseholderOnTheLeft(
H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
120 t.noalias() = mat * v;
121 v = precond.solve(t);
125 for (Index i = 0; i < k; ++i) {
126 v.tail(m - i).applyHouseholderOnTheLeft(
H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
129 if (v.tail(m - k).norm() != 0.0)
134 Ref<VectorType> Hk_tail =
H.col(k).tail(m - k - 1);
135 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
138 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
144 for (Index i = 0; i < k - 1; ++i)
147 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
151 if (k<m && v(k) != (Scalar) 0)
154 G[k - 1].makeGivens(v(k - 1), v(k));
157 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
158 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
162 H.col(k-1).head(k) = v.head(k);
164 tol_error = abs(w(k)) / r0Norm;
165 std::cout <<
"iteration: " << iters <<
", error: " << tol_error << std::endl;
166 bool stop = (k==m || tol_error < tol || iters == maxIters);
168 if (stop || k == restart)
171 Ref<VectorType> y = w.head(k);
172 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
176 for (Index i = k - 1; i >= 0; --i)
180 x_new.tail(m - i).applyHouseholderOnTheLeft(
H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
194 p0.noalias() = rhs - mat*x;
195 r0 = precond.solve(p0);
203 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
215 template<
typename _MatrixType,
216 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
221 template<
typename _MatrixType,
typename _Preconditioner>
222 struct traits<
GMRES<_MatrixType,_Preconditioner> >
224 typedef _MatrixType MatrixType;
225 typedef _Preconditioner Preconditioner;
264 template<
typename _MatrixType,
typename _Preconditioner>
265 class GMRES :
public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
267 typedef IterativeSolverBase<GMRES> Base;
270 using Base::m_iterations;
272 using Base::m_isInitialized;
278 using Base::_solve_impl;
279 typedef _MatrixType MatrixType;
280 typedef typename MatrixType::Scalar Scalar;
281 typedef typename MatrixType::RealScalar RealScalar;
282 typedef _Preconditioner Preconditioner;
299 template<
typename MatrixDerived>
300 explicit GMRES(
const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
314 template<
typename Rhs,
typename Dest>
315 void _solve_with_guess_impl(
const Rhs& b, Dest& x)
const
318 for(Index j=0; j<b.cols(); ++j)
320 m_iterations = Base::maxIterations();
321 m_error = Base::m_tolerance;
323 typename Dest::ColXpr xj(x,j);
324 if(!internal::gmres(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
327 m_info = failed ? NumericalIssue
328 : m_error <= Base::m_tolerance ? Success
330 m_isInitialized =
true;
334 template<
typename Rhs,
typename Dest>
335 void _solve_impl(
const Rhs& b, MatrixBase<Dest> &x)
const
338 if(x.squaredNorm() == 0)
return;
339 _solve_with_guess_impl(b,x.derived());
348 #endif // EIGEN_GMRES_H