13#include "../../../../Eigen/Eigenvalues"
16#include "./InternalHeaderCheck.h"
20template <
typename MatrixType_,
typename Preconditioner_ = DiagonalPreconditioner<
typename MatrixType_::Scalar> >
25template <
typename MatrixType_,
typename Preconditioner_>
26struct traits<DGMRES<MatrixType_, Preconditioner_> > {
27 typedef MatrixType_ MatrixType;
28 typedef Preconditioner_ Preconditioner;
39template <
typename VectorType,
typename IndexType>
40void sortWithPermutation(VectorType& vec, IndexType& perm,
typename IndexType::Scalar& ncut) {
41 eigen_assert(vec.size() == perm.size());
43 for (
Index k = 0; k < ncut; k++) {
45 for (
Index j = 0; j < vec.size() - 1; j++) {
46 if (vec(perm(j)) < vec(perm(j + 1))) {
47 std::swap(perm(j), perm(j + 1));
97template <
typename MatrixType_,
typename Preconditioner_>
102 using Base::m_isInitialized;
103 using Base::m_iterations;
104 using Base::m_tolerance;
108 using Base::_solve_impl;
109 using Base::_solve_with_guess_impl;
110 typedef MatrixType_ MatrixType;
111 typedef typename MatrixType::Scalar Scalar;
112 typedef typename MatrixType::StorageIndex StorageIndex;
113 typedef typename MatrixType::RealScalar RealScalar;
114 typedef internal::make_complex_t<Scalar> ComplexScalar;
115 typedef Preconditioner_ Preconditioner;
124 : Base(), m_restart(30), m_neig(0), m_r(0), m_maxNeig(5), m_isDeflAllocated(false), m_isDeflInitialized(false) {}
136 template <
typename MatrixDerived>
143 m_isDeflAllocated(false),
144 m_isDeflInitialized(false) {}
149 template <
typename Rhs,
typename Dest>
150 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
151 EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime == 1 || Dest::ColsAtCompileTime == 1,
152 YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
155 m_error = Base::m_tolerance;
157 dgmres(matrix(), b, x, Base::m_preconditioner);
175 if (neig + 1 > m_maxNeig) m_maxNeig = neig + 1;
190 template <
typename Rhs,
typename Dest>
191 void dgmres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond)
const;
193 template <
typename Dest>
194 Index dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
195 const RealScalar& normRhs,
Index& nbIts)
const;
197 Index dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it,
198 StorageIndex& neig)
const;
200 template <
typename RhsType,
typename DestType>
201 Index dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
205 void dgmresInitDeflation(
Index& rows)
const;
206 mutable DenseMatrix m_V;
207 mutable DenseMatrix m_H;
208 mutable DenseMatrix m_Hes;
209 mutable Index m_restart;
210 mutable DenseMatrix m_U;
211 mutable DenseMatrix m_MU;
212 mutable DenseMatrix m_T;
214 mutable StorageIndex m_neig;
216 mutable Index m_maxNeig;
217 mutable RealScalar m_lambdaN;
218 mutable bool m_isDeflAllocated;
219 mutable bool m_isDeflInitialized;
222 mutable RealScalar m_smv;
223 mutable bool m_force;
231template <
typename MatrixType_,
typename Preconditioner_>
232template <
typename Rhs,
typename Dest>
234 const Preconditioner& precond)
const {
235 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
237 RealScalar normRhs = rhs.norm();
238 if (normRhs <= considerAsZero) {
245 m_isDeflInitialized =
false;
246 Index n = mat.rows();
249 m_H.resize(m_restart + 1, m_restart);
250 m_Hes.resize(m_restart, m_restart);
251 m_V.resize(n, m_restart + 1);
253 if (x.squaredNorm() == 0) x = precond.solve(rhs);
255 RealScalar beta = r0.norm();
257 m_error = beta / normRhs;
258 if (m_error < m_tolerance)
265 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
285template <
typename MatrixType_,
typename Preconditioner_>
286template <
typename Dest>
288 DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
289 Index& nbIts)
const {
291 DenseVector g(m_restart + 1);
294 m_V.col(0) = r0 / beta;
296 std::vector<JacobiRotation<Scalar> > gr(m_restart);
298 Index n = mat.rows();
299 DenseVector tv1(n), tv2(n);
300 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations) {
302 if (m_isDeflInitialized) {
303 dgmresApplyDeflation(m_V.col(it), tv1);
304 tv2 = precond.solve(tv1);
306 tv2 = precond.solve(m_V.col(it));
312 for (
Index i = 0; i <= it; ++i) {
313 coef = tv1.dot(m_V.col(i));
314 tv1 = tv1 - coef * m_V.col(i);
320 m_V.col(it + 1) = tv1 / coef;
321 m_H(it + 1, it) = coef;
327 for (
Index i = 1; i <= it; ++i) {
328 m_H.col(it).applyOnTheLeft(i - 1, i, gr[i - 1].adjoint());
331 gr[it].makeGivens(m_H(it, it), m_H(it + 1, it));
333 m_H.col(it).applyOnTheLeft(it, it + 1, gr[it].adjoint());
334 g.applyOnTheLeft(it, it + 1, gr[it].adjoint());
336 beta = std::abs(g(it + 1));
337 m_error = beta / normRhs;
342 if (m_error < m_tolerance) {
352 DenseVector nrs(m_restart);
353 nrs = m_H.topLeftCorner(it, it).template triangularView<Upper>().solve(g.head(it));
356 if (m_isDeflInitialized) {
357 tv1 = m_V.leftCols(it) * nrs;
358 dgmresApplyDeflation(tv1, tv2);
359 x = x + precond.solve(tv2);
361 x = x + precond.solve(m_V.leftCols(it) * nrs);
364 if (nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r + m_neig) < m_maxNeig)
365 dgmresComputeDeflationData(mat, precond, it, m_neig);
369template <
typename MatrixType_,
typename Preconditioner_>
370void DGMRES<MatrixType_, Preconditioner_>::dgmresInitDeflation(
Index& rows)
const {
371 m_U.resize(rows, m_maxNeig);
372 m_MU.resize(rows, m_maxNeig);
373 m_T.resize(m_maxNeig, m_maxNeig);
375 m_isDeflAllocated =
true;
378template <
typename MatrixType_,
typename Preconditioner_>
379inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
380 const ComplexSchur<DenseMatrix>& schurofH)
const {
381 return schurofH.matrixT().diagonal();
384template <
typename MatrixType_,
typename Preconditioner_>
385inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
387 const DenseMatrix& T = schurofH.matrixT();
389 ComplexVector eig(it);
392 if (T(j + 1, j) ==
Scalar(0)) {
393 eig(j) = ComplexScalar(T(j, j), RealScalar(0));
396 eig(j) = ComplexScalar(T(j, j), T(j + 1, j));
397 eig(j + 1) = ComplexScalar(T(j, j + 1), T(j + 1, j + 1));
401 if (j < it - 1) eig(j) = ComplexScalar(T(j, j), RealScalar(0));
405template <
typename MatrixType_,
typename Preconditioner_>
406Index DGMRES<MatrixType_, Preconditioner_>::dgmresComputeDeflationData(
const MatrixType& mat,
407 const Preconditioner& precond,
const Index& it,
408 StorageIndex& neig)
const {
411 bool computeU =
true;
412 DenseMatrix matrixQ(it, it);
413 matrixQ.setIdentity();
414 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it, it), matrixQ, computeU);
416 ComplexVector eig(it);
418 eig = this->schurValues(schurofH);
421 DenseRealVector modulEig(it);
422 for (
Index j = 0; j < it; ++j) modulEig(j) = std::abs(eig(j));
423 perm.setLinSpaced(it, 0, internal::convert_index<StorageIndex>(it - 1));
424 internal::sortWithPermutation(modulEig, perm, neig);
427 m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
431 while (nbrEig < neig) {
432 if (eig(perm(it - nbrEig - 1)).
imag() == RealScalar(0))
438 DenseMatrix Sr(it, nbrEig);
440 for (
Index j = 0; j < nbrEig; j++) {
441 Sr.col(j) = schurofH.matrixU().col(perm(it - j - 1));
446 X = m_V.leftCols(it) * Sr;
449 for (
Index j = 0; j < nbrEig; j++)
450 for (
Index k = 0; k < m_r; k++) X.col(j) = X.col(j) - (m_U.col(k).dot(X.col(j))) * m_U.col(k);
454 Index m = m_V.rows();
455 if (!m_isDeflAllocated) dgmresInitDeflation(m);
456 DenseMatrix MX(m, nbrEig);
458 for (
Index j = 0; j < nbrEig; j++) {
459 tv1 = mat * X.col(j);
460 MX.col(j) = precond.solve(tv1);
464 m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
466 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
467 m_T.block(m_r, 0, nbrEig, m_r) = X.transpose() * m_MU.leftCols(m_r);
471 for (
Index j = 0; j < nbrEig; j++) m_U.col(m_r + j) = X.col(j);
472 for (
Index j = 0; j < nbrEig; j++) m_MU.col(m_r + j) = MX.col(j);
477 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
480 m_isDeflInitialized =
true;
483template <
typename MatrixType_,
typename Preconditioner_>
484template <
typename RhsType,
typename DestType>
485Index DGMRES<MatrixType_, Preconditioner_>::dgmresApplyDeflation(
const RhsType& x, DestType& y)
const {
486 DenseVector x1 = m_U.leftCols(m_r).transpose() * x;
487 y = x + m_U.leftCols(m_r) * (m_lambdaN * m_luT.solve(x1) - x1);
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition DGMRES.h:98
Index dgmresCycle(const MatrixType &mat, const Preconditioner &precond, Dest &x, DenseVector &r0, RealScalar &beta, const RealScalar &normRhs, Index &nbIts) const
Perform one restart cycle of DGMRES.
Definition DGMRES.h:287
Index restart()
Definition DGMRES.h:163
DGMRES(const EigenBase< MatrixDerived > &A)
Definition DGMRES.h:137
void setMaxEigenv(const Index maxNeig)
Definition DGMRES.h:186
DGMRES()
Definition DGMRES.h:123
void dgmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond) const
Perform several cycles of restarted GMRES with modified Gram Schmidt,.
Definition DGMRES.h:233
void set_restart(const Index restart)
Definition DGMRES.h:168
Index deflSize()
Definition DGMRES.h:181
void setEigenv(const Index neig)
Definition DGMRES.h:173
Index maxIterations() const
Derived & setZero(Index rows, Index cols)
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_imag_op< typename Derived::Scalar >, const Derived > imag(const Eigen::ArrayBase< Derived > &x)