14 #include "./InternalHeaderCheck.h"
57 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
58 bool gmres(
const MatrixType & mat,
const Rhs & rhs, Dest & x,
const Preconditioner & precond,
59 Index &iters,
const Index &restart,
typename Dest::RealScalar & tol_error) {
64 typedef typename Dest::RealScalar RealScalar;
65 typedef typename Dest::Scalar Scalar;
67 typedef Matrix < Scalar, Dynamic, Dynamic, ColMajor> FMatrixType;
69 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
71 if(rhs.norm() <= considerAsZero)
78 RealScalar tol = tol_error;
79 const Index maxIters = iters;
82 const Index m = mat.rows();
85 VectorType p0 = rhs - mat*x;
86 VectorType r0 = precond.solve(p0);
88 const RealScalar r0Norm = r0.norm();
98 FMatrixType H = FMatrixType::Zero(m, restart + 1);
99 VectorType w = VectorType::Zero(restart + 1);
100 VectorType tau = VectorType::Zero(restart + 1);
103 std::vector < JacobiRotation < Scalar > > G(restart);
106 VectorType t(m), v(m), workspace(m), x_new(m);
109 Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
111 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
114 for (
Index k = 1; k <= restart; ++k)
118 v = VectorType::Unit(m, k - 1);
122 for (
Index i = k - 1; i >= 0; --i) {
123 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
127 t.noalias() = mat * v;
128 v = precond.solve(t);
132 for (
Index i = 0; i < k; ++i) {
133 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
136 if (v.tail(m - k).norm() != 0.0)
141 Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
142 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
145 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
151 for (
Index i = 0; i < k - 1; ++i)
154 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
158 if (k<m && v(k) != (Scalar) 0)
161 G[k - 1].makeGivens(v(k - 1), v(k));
164 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
165 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
169 H.col(k-1).head(k) = v.head(k);
171 tol_error =
abs(w(k)) / r0Norm;
172 bool stop = (k==m || tol_error < tol || iters == maxIters);
174 if (stop || k == restart)
177 Ref<VectorType> y = w.head(k);
178 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
182 for (
Index i = k - 1; i >= 0; --i)
186 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
200 p0.noalias() = rhs - mat*x;
201 r0 = precond.solve(p0);
209 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
221 template<
typename MatrixType_,
222 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
227 template<
typename MatrixType_,
typename Preconditioner_>
228 struct traits<GMRES<MatrixType_,Preconditioner_> >
230 typedef MatrixType_ MatrixType;
231 typedef Preconditioner_ Preconditioner;
270 template<
typename MatrixType_,
typename Preconditioner_>
276 using Base::m_iterations;
278 using Base::m_isInitialized;
284 using Base::_solve_impl;
285 typedef MatrixType_ MatrixType;
286 typedef typename MatrixType::Scalar Scalar;
287 typedef typename MatrixType::RealScalar RealScalar;
288 typedef Preconditioner_ Preconditioner;
305 template<
typename MatrixDerived>
320 template<
typename Rhs,
typename Dest>
321 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
324 m_error = Base::m_tolerance;
325 bool ret = internal::gmres(matrix(), b, x, Base::m_preconditioner, m_iterations, m_restart, m_error);
327 : m_error <= Base::m_tolerance ?
Success
A GMRES solver for sparse square problems.
Definition: GMRES.h:272
GMRES(const EigenBase< MatrixDerived > &A)
Definition: GMRES.h:306
Index get_restart()
Definition: GMRES.h:312
GMRES()
Definition: GMRES.h:293
void set_restart(const Index restart)
Definition: GMRES.h:317
Index maxIterations() const
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)