10 #ifndef EIGEN_DGMRES_H
11 #define EIGEN_DGMRES_H
13 #include "../../../../Eigen/Eigenvalues"
15 #include "./InternalHeaderCheck.h"
19 template<
typename MatrixType_,
20 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
25 template<
typename MatrixType_,
typename Preconditioner_>
26 struct traits<DGMRES<MatrixType_,Preconditioner_> >
28 typedef MatrixType_ MatrixType;
29 typedef Preconditioner_ Preconditioner;
40 template <
typename VectorType,
typename IndexType>
41 void sortWithPermutation (VectorType& vec, IndexType& perm,
typename IndexType::Scalar& ncut)
43 eigen_assert(vec.size() == perm.size());
45 for (
Index k = 0; k < ncut; k++)
48 for (
Index j = 0; j < vec.size()-1; j++)
50 if ( vec(perm(j)) < vec(perm(j+1)) )
52 std::swap(perm(j),perm(j+1));
102 template<
typename MatrixType_,
typename Preconditioner_>
108 using Base::m_iterations;
110 using Base::m_isInitialized;
111 using Base::m_tolerance;
113 using Base::_solve_impl;
114 using Base::_solve_with_guess_impl;
115 typedef MatrixType_ MatrixType;
116 typedef typename MatrixType::Scalar Scalar;
117 typedef typename MatrixType::StorageIndex StorageIndex;
118 typedef typename MatrixType::RealScalar RealScalar;
119 typedef Preconditioner_ Preconditioner;
128 DGMRES() :
Base(),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
140 template<
typename MatrixDerived>
146 template<
typename Rhs,
typename Dest>
147 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
149 EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
152 m_error = Base::m_tolerance;
154 dgmres(matrix(), b, x, Base::m_preconditioner);
173 if (neig+1 > m_maxNeig) m_maxNeig = neig+1;
188 template<
typename Rhs,
typename Dest>
189 void dgmres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond)
const;
191 template<
typename Dest>
192 Index dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
Index& nbIts)
const;
194 Index dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, StorageIndex& neig)
const;
196 template<
typename RhsType,
typename DestType>
197 Index dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
201 void dgmresInitDeflation(
Index& rows)
const;
202 mutable DenseMatrix m_V;
203 mutable DenseMatrix m_H;
204 mutable DenseMatrix m_Hes;
205 mutable Index m_restart;
206 mutable DenseMatrix m_U;
207 mutable DenseMatrix m_MU;
208 mutable DenseMatrix m_T;
210 mutable StorageIndex m_neig;
212 mutable Index m_maxNeig;
213 mutable RealScalar m_lambdaN;
214 mutable bool m_isDeflAllocated;
215 mutable bool m_isDeflInitialized;
218 mutable RealScalar m_smv;
219 mutable bool m_force;
228 template<
typename MatrixType_,
typename Preconditioner_>
229 template<
typename Rhs,
typename Dest>
231 const Preconditioner& precond)
const
233 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
235 RealScalar normRhs = rhs.norm();
236 if(normRhs <= considerAsZero)
244 m_isDeflInitialized =
false;
245 Index n = mat.rows();
248 m_H.resize(m_restart+1, m_restart);
249 m_Hes.resize(m_restart, m_restart);
250 m_V.resize(n,m_restart+1);
252 if(x.squaredNorm()==0)
253 x = precond.solve(rhs);
255 RealScalar beta = r0.norm();
257 m_error = beta/normRhs;
258 if(m_error < m_tolerance)
266 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
286 template<
typename MatrixType_,
typename Preconditioner_>
287 template<
typename Dest>
294 m_V.col(0) = r0/beta;
296 std::vector<JacobiRotation<Scalar> >gr(m_restart);
298 Index n = mat.rows();
300 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
303 if (m_isDeflInitialized )
305 dgmresApplyDeflation(m_V.col(it), tv1);
306 tv2 = precond.solve(tv1);
310 tv2 = precond.solve(m_V.col(it));
316 for (
Index i = 0; i <= it; ++i)
318 coef = tv1.dot(m_V.col(i));
319 tv1 = tv1 - coef * m_V.col(i);
325 m_V.col(it+1) = tv1/coef;
326 m_H(it+1, it) = coef;
332 for (
Index i = 1; i <= it; ++i)
334 m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
337 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
339 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
340 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
342 beta = std::abs(g(it+1));
343 m_error = beta/normRhs;
347 if (m_error < m_tolerance)
359 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
362 if (m_isDeflInitialized)
364 tv1 = m_V.leftCols(it) * nrs;
365 dgmresApplyDeflation(tv1, tv2);
366 x = x + precond.solve(tv2);
369 x = x + precond.solve(m_V.leftCols(it) * nrs);
372 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
373 dgmresComputeDeflationData(mat, precond, it, m_neig);
379 template<
typename MatrixType_,
typename Preconditioner_>
382 m_U.resize(rows, m_maxNeig);
383 m_MU.resize(rows, m_maxNeig);
384 m_T.resize(m_maxNeig, m_maxNeig);
386 m_isDeflAllocated =
true;
389 template<
typename MatrixType_,
typename Preconditioner_>
390 inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
const ComplexSchur<DenseMatrix>& schurofH)
const
392 return schurofH.matrixT().diagonal();
395 template<
typename MatrixType_,
typename Preconditioner_>
396 inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
const RealSchur<DenseMatrix>& schurofH)
const
398 const DenseMatrix& T = schurofH.matrixT();
400 ComplexVector eig(it);
404 if (T(j+1,j) ==Scalar(0))
406 eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
411 eig(j) = std::complex<RealScalar>(T(j,j),T(j+1,j));
412 eig(j+1) = std::complex<RealScalar>(T(j,j+1),T(j+1,j+1));
416 if (j < it-1) eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
420 template<
typename MatrixType_,
typename Preconditioner_>
421 Index DGMRES<MatrixType_, Preconditioner_>::dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, StorageIndex& neig)
const
424 std::conditional_t<NumTraits<Scalar>::IsComplex, ComplexSchur<DenseMatrix>, RealSchur<DenseMatrix> > schurofH;
425 bool computeU =
true;
426 DenseMatrix matrixQ(it,it);
427 matrixQ.setIdentity();
428 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
430 ComplexVector eig(it);
432 eig = this->schurValues(schurofH);
435 DenseRealVector modulEig(it);
436 for (
Index j=0; j<it; ++j) modulEig(j) = std::abs(eig(j));
437 perm.setLinSpaced(it,0,internal::convert_index<StorageIndex>(it-1));
438 internal::sortWithPermutation(modulEig, perm, neig);
442 m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
446 while (nbrEig < neig)
448 if(eig(perm(it-nbrEig-1)).imag() == RealScalar(0)) nbrEig++;
452 DenseMatrix Sr(it, nbrEig);
454 for (
Index j = 0; j < nbrEig; j++)
456 Sr.col(j) = schurofH.matrixU().col(perm(it-j-1));
461 X = m_V.leftCols(it) * Sr;
465 for (
Index j = 0; j < nbrEig; j++)
466 for (
Index k =0; k < m_r; k++)
467 X.col(j) = X.col(j) - (m_U.col(k).dot(X.col(j)))*m_U.col(k);
471 Index m = m_V.rows();
472 if (!m_isDeflAllocated)
473 dgmresInitDeflation(m);
474 DenseMatrix MX(m, nbrEig);
476 for (
Index j = 0; j < nbrEig; j++)
478 tv1 = mat * X.col(j);
479 MX.col(j) = precond.solve(tv1);
483 m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
486 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
487 m_T.block(m_r, 0, nbrEig, m_r) = X.transpose() * m_MU.leftCols(m_r);
491 for (
Index j = 0; j < nbrEig; j++) m_U.col(m_r+j) = X.col(j);
492 for (
Index j = 0; j < nbrEig; j++) m_MU.col(m_r+j) = MX.col(j);
497 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
500 m_isDeflInitialized =
true;
503 template<
typename MatrixType_,
typename Preconditioner_>
504 template<
typename RhsType,
typename DestType>
505 Index DGMRES<MatrixType_, Preconditioner_>::dgmresApplyDeflation(
const RhsType &x, DestType &y)
const
507 DenseVector x1 = m_U.leftCols(m_r).transpose() * x;
508 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:104
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:288
Index restart()
Definition: DGMRES.h:160
DGMRES(const EigenBase< MatrixDerived > &A)
Definition: DGMRES.h:141
void setMaxEigenv(const Index maxNeig)
Definition: DGMRES.h:184
DGMRES()
Definition: DGMRES.h:128
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:230
void set_restart(const Index restart)
Definition: DGMRES.h:165
Index deflSize()
Definition: DGMRES.h:179
void setEigenv(const Index neig)
Definition: DGMRES.h:170
Index maxIterations() const
Derived & setZero(Index size)
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index