10 #ifndef EIGEN_DGMRES_H 11 #define EIGEN_DGMRES_H 13 #include "../../../../Eigen/Eigenvalues" 17 template<
typename _MatrixType,
18 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
23 template<
typename _MatrixType,
typename _Preconditioner>
24 struct traits<
DGMRES<_MatrixType,_Preconditioner> >
26 typedef _MatrixType MatrixType;
27 typedef _Preconditioner Preconditioner;
38 template <
typename VectorType,
typename IndexType>
39 void sortWithPermutation (VectorType& vec, IndexType& perm,
typename IndexType::Scalar& ncut)
41 eigen_assert(vec.size() == perm.size());
43 for (
Index k = 0; k < ncut; k++)
46 for (
Index j = 0; j < vec.size()-1; j++)
48 if ( vec(perm(j)) < vec(perm(j+1)) )
50 std::swap(perm(j),perm(j+1));
100 template<
typename _MatrixType,
typename _Preconditioner>
106 using Base::m_iterations;
108 using Base::m_isInitialized;
109 using Base::m_tolerance;
111 using Base::_solve_impl;
112 using Base::_solve_with_guess_impl;
113 typedef _MatrixType MatrixType;
114 typedef typename MatrixType::Scalar Scalar;
115 typedef typename MatrixType::StorageIndex StorageIndex;
116 typedef typename MatrixType::RealScalar RealScalar;
117 typedef _Preconditioner Preconditioner;
126 DGMRES() : Base(),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
138 template<
typename MatrixDerived>
139 explicit DGMRES(
const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
144 template<
typename Rhs,
typename Dest>
145 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const 147 EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
149 m_iterations = Base::maxIterations();
150 m_error = Base::m_tolerance;
152 dgmres(matrix(), b, x, Base::m_preconditioner);
171 if (neig+1 > m_maxNeig) m_maxNeig = neig+1;
186 template<
typename Rhs,
typename Dest>
187 void dgmres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond)
const;
189 template<
typename Dest>
190 Index dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
Index& nbIts)
const;
192 Index dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, StorageIndex& neig)
const;
194 template<
typename RhsType,
typename DestType>
195 Index dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
199 void dgmresInitDeflation(
Index& rows)
const;
200 mutable DenseMatrix m_V;
201 mutable DenseMatrix m_H;
202 mutable DenseMatrix m_Hes;
203 mutable Index m_restart;
204 mutable DenseMatrix m_U;
205 mutable DenseMatrix m_MU;
206 mutable DenseMatrix m_T;
208 mutable StorageIndex m_neig;
210 mutable Index m_maxNeig;
211 mutable RealScalar m_lambdaN;
212 mutable bool m_isDeflAllocated;
213 mutable bool m_isDeflInitialized;
216 mutable RealScalar m_smv;
217 mutable bool m_force;
226 template<
typename _MatrixType,
typename _Preconditioner>
227 template<
typename Rhs,
typename Dest>
229 const Preconditioner& precond)
const 231 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
233 RealScalar normRhs = rhs.norm();
234 if(normRhs <= considerAsZero)
242 m_isDeflInitialized =
false;
243 Index n = mat.rows();
246 m_H.resize(m_restart+1, m_restart);
247 m_Hes.resize(m_restart, m_restart);
248 m_V.resize(n,m_restart+1);
250 if(x.squaredNorm()==0)
251 x = precond.solve(rhs);
253 RealScalar beta = r0.norm();
255 m_error = beta/normRhs;
256 if(m_error < m_tolerance)
264 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
284 template<
typename _MatrixType,
typename _Preconditioner>
285 template<
typename Dest>
292 m_V.col(0) = r0/beta;
294 std::vector<JacobiRotation<Scalar> >gr(m_restart);
296 Index n = mat.rows();
298 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
301 if (m_isDeflInitialized )
303 dgmresApplyDeflation(m_V.col(it), tv1);
304 tv2 = precond.solve(tv1);
308 tv2 = precond.solve(m_V.col(it));
314 for (
Index i = 0; i <= it; ++i)
316 coef = tv1.dot(m_V.col(i));
317 tv1 = tv1 - coef * m_V.col(i);
323 m_V.col(it+1) = tv1/coef;
324 m_H(it+1, it) = coef;
330 for (
Index i = 1; i <= it; ++i)
332 m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
335 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
337 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
338 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
340 beta = std::abs(g(it+1));
341 m_error = beta/normRhs;
345 if (m_error < m_tolerance)
357 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
360 if (m_isDeflInitialized)
362 tv1 = m_V.leftCols(it) * nrs;
363 dgmresApplyDeflation(tv1, tv2);
364 x = x + precond.solve(tv2);
367 x = x + precond.solve(m_V.leftCols(it) * nrs);
370 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
371 dgmresComputeDeflationData(mat, precond, it, m_neig);
377 template<
typename _MatrixType,
typename _Preconditioner>
380 m_U.resize(rows, m_maxNeig);
381 m_MU.resize(rows, m_maxNeig);
382 m_T.resize(m_maxNeig, m_maxNeig);
384 m_isDeflAllocated =
true;
387 template<
typename _MatrixType,
typename _Preconditioner>
390 return schurofH.
matrixT().diagonal();
393 template<
typename _MatrixType,
typename _Preconditioner>
402 if (T(j+1,j) ==Scalar(0))
404 eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
409 eig(j) = std::complex<RealScalar>(T(j,j),T(j+1,j));
410 eig(j+1) = std::complex<RealScalar>(T(j,j+1),T(j+1,j+1));
414 if (j < it-1) eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
418 template<
typename _MatrixType,
typename _Preconditioner>
423 bool computeU =
true;
425 matrixQ.setIdentity();
426 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
430 eig = this->schurValues(schurofH);
434 for (
Index j=0; j<it; ++j) modulEig(j) = std::abs(eig(j));
435 perm.setLinSpaced(it,0,internal::convert_index<StorageIndex>(it-1));
436 internal::sortWithPermutation(modulEig, perm, neig);
440 m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
444 while (nbrEig < neig)
446 if(eig(perm(it-nbrEig-1)).
imag() == RealScalar(0)) nbrEig++;
452 for (
Index j = 0; j < nbrEig; j++)
454 Sr.col(j) = schurofH.matrixU().col(perm(it-j-1));
459 X = m_V.leftCols(it) * Sr;
463 for (
Index j = 0; j < nbrEig; j++)
464 for (
Index k =0; k < m_r; k++)
465 X.col(j) = X.col(j) - (m_U.col(k).dot(X.col(j)))*m_U.col(k);
469 Index m = m_V.rows();
470 if (!m_isDeflAllocated)
471 dgmresInitDeflation(m);
474 for (
Index j = 0; j < nbrEig; j++)
476 tv1 = mat * X.col(j);
477 MX.col(j) = precond.solve(tv1);
481 m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
484 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
485 m_T.block(m_r, 0, nbrEig, m_r) = X.transpose() * m_MU.leftCols(m_r);
489 for (
Index j = 0; j < nbrEig; j++) m_U.col(m_r+j) = X.col(j);
490 for (
Index j = 0; j < nbrEig; j++) m_MU.col(m_r+j) = MX.col(j);
495 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
498 m_isDeflInitialized =
true;
501 template<
typename _MatrixType,
typename _Preconditioner>
502 template<
typename RhsType,
typename DestType>
505 DenseVector x1 = m_U.leftCols(m_r).transpose() * x;
506 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:19
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:228
void setEigenv(const Index neig)
Definition: DGMRES.h:168
Index deflSize()
Definition: DGMRES.h:177
Namespace containing all symbols from the Eigen library.
Derived & setZero(Index size)
void set_restart(const Index restart)
Definition: DGMRES.h:163
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const ImagReturnType imag() const
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:286
Index restart()
Definition: DGMRES.h:158
const MatrixType & matrixT() const
DGMRES()
Definition: DGMRES.h:126
void setMaxEigenv(const Index maxNeig)
Definition: DGMRES.h:182
const ComplexMatrixType & matrixT() const
DGMRES(const EigenBase< MatrixDerived > &A)
Definition: DGMRES.h:139