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());
42 typedef typename IndexType::Scalar Index;
43 typedef typename VectorType::Scalar Scalar;
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));
100 template<
typename _MatrixType,
typename _Preconditioner>
101 class DGMRES :
public IterativeSolverBase<DGMRES<_MatrixType,_Preconditioner> >
103 typedef IterativeSolverBase<DGMRES> Base;
104 using Base::mp_matrix;
106 using Base::m_iterations;
108 using Base::m_isInitialized;
109 using Base::m_tolerance;
111 typedef _MatrixType MatrixType;
112 typedef typename MatrixType::Scalar Scalar;
113 typedef typename MatrixType::Index Index;
114 typedef typename MatrixType::RealScalar RealScalar;
115 typedef _Preconditioner Preconditioner;
116 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
117 typedef Matrix<RealScalar,Dynamic,Dynamic> DenseRealMatrix;
118 typedef Matrix<Scalar,Dynamic,1> DenseVector;
119 typedef Matrix<RealScalar,Dynamic,1> DenseRealVector;
120 typedef Matrix<std::complex<RealScalar>, Dynamic, 1> ComplexVector;
124 DGMRES() : Base(),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
136 template<
typename MatrixDerived>
137 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) {}
146 template<
typename Rhs,
typename Guess>
147 inline const internal::solve_retval_with_guess<DGMRES, Rhs, Guess>
150 eigen_assert(m_isInitialized &&
"DGMRES is not initialized.");
151 eigen_assert(Base::rows()==b.rows()
152 &&
"DGMRES::solve(): invalid number of rows of the right hand side matrix b");
153 return internal::solve_retval_with_guess
154 <
DGMRES, Rhs, Guess>(*
this, b.derived(), x0);
158 template<
typename Rhs,
typename Dest>
159 void _solveWithGuess(
const Rhs& b, Dest& x)
const
162 for(
int j=0; j<b.cols(); ++j)
164 m_iterations = Base::maxIterations();
165 m_error = Base::m_tolerance;
167 typename Dest::ColXpr xj(x,j);
168 dgmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner);
170 m_info = failed ? NumericalIssue
171 : m_error <= Base::m_tolerance ? Success
173 m_isInitialized =
true;
177 template<
typename Rhs,
typename Dest>
178 void _solve(
const Rhs& b, Dest& x)
const
181 _solveWithGuess(b,x);
199 if (neig+1 > m_maxNeig) m_maxNeig = neig+1;
214 template<
typename Rhs,
typename Dest>
215 void dgmres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond)
const;
217 template<
typename Dest>
218 int dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
int& nbIts)
const;
220 int dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, Index& neig)
const;
222 template<
typename RhsType,
typename DestType>
223 int dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
224 ComplexVector schurValues(
const ComplexSchur<DenseMatrix>& schurofH)
const;
225 ComplexVector schurValues(
const RealSchur<DenseMatrix>& schurofH)
const;
227 void dgmresInitDeflation(Index& rows)
const;
228 mutable DenseMatrix m_V;
229 mutable DenseMatrix m_H;
230 mutable DenseMatrix m_Hes;
231 mutable Index m_restart;
232 mutable DenseMatrix m_U;
233 mutable DenseMatrix m_MU;
234 mutable DenseMatrix m_T;
235 mutable PartialPivLU<DenseMatrix> m_luT;
238 mutable int m_maxNeig;
239 mutable RealScalar m_lambdaN;
240 mutable bool m_isDeflAllocated;
241 mutable bool m_isDeflInitialized;
244 mutable RealScalar m_smv;
245 mutable bool m_force;
254 template<
typename _MatrixType,
typename _Preconditioner>
255 template<
typename Rhs,
typename Dest>
257 const Preconditioner& precond)
const
263 m_H.resize(m_restart+1, m_restart);
264 m_Hes.resize(m_restart, m_restart);
265 m_V.resize(n,m_restart+1);
267 x = precond.solve(x);
269 RealScalar beta = r0.norm();
270 RealScalar normRhs = rhs.norm();
271 m_error = beta/normRhs;
272 if(m_error < m_tolerance)
275 m_info = NoConvergence;
278 while (nbIts < m_iterations && m_info == NoConvergence)
280 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
283 if (nbIts < m_iterations && m_info == NoConvergence)
298 template<
typename _MatrixType,
typename _Preconditioner>
299 template<
typename Dest>
303 DenseVector g(m_restart+1);
306 m_V.col(0) = r0/beta;
307 m_info = NoConvergence;
308 std::vector<JacobiRotation<Scalar> >gr(m_restart);
311 DenseVector tv1(n), tv2(n);
312 while (m_info == NoConvergence && it < m_restart && nbIts < m_iterations)
315 if (m_isDeflInitialized )
317 dgmresApplyDeflation(m_V.col(it), tv1);
318 tv2 = precond.solve(tv1);
322 tv2 = precond.solve(m_V.col(it));
328 for (
int i = 0; i <= it; ++i)
330 coef = tv1.dot(m_V.col(i));
331 tv1 = tv1 - coef * m_V.col(i);
337 m_V.col(it+1) = tv1/coef;
338 m_H(it+1, it) = coef;
344 for (
int i = 1; i <= it; ++i)
346 m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
349 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
351 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
352 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
354 beta = std::abs(g(it+1));
355 m_error = beta/normRhs;
356 std::cerr << nbIts <<
" Relative Residual Norm " << m_error << std::endl;
359 if (m_error < m_tolerance)
370 DenseVector nrs(m_restart);
371 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
374 if (m_isDeflInitialized)
376 tv1 = m_V.leftCols(it) * nrs;
377 dgmresApplyDeflation(tv1, tv2);
378 x = x + precond.solve(tv2);
381 x = x + precond.solve(m_V.leftCols(it) * nrs);
384 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
385 dgmresComputeDeflationData(mat, precond, it, m_neig);
391 template<
typename _MatrixType,
typename _Preconditioner>
394 m_U.resize(rows, m_maxNeig);
395 m_MU.resize(rows, m_maxNeig);
396 m_T.resize(m_maxNeig, m_maxNeig);
398 m_isDeflAllocated =
true;
401 template<
typename _MatrixType,
typename _Preconditioner>
402 inline typename DGMRES<_MatrixType, _Preconditioner>::ComplexVector DGMRES<_MatrixType, _Preconditioner>::schurValues(
const ComplexSchur<DenseMatrix>& schurofH)
const
404 return schurofH.matrixT().diagonal();
407 template<
typename _MatrixType,
typename _Preconditioner>
408 inline typename DGMRES<_MatrixType, _Preconditioner>::ComplexVector DGMRES<_MatrixType, _Preconditioner>::schurValues(
const RealSchur<DenseMatrix>& schurofH)
const
410 typedef typename MatrixType::Index Index;
411 const DenseMatrix& T = schurofH.matrixT();
413 ComplexVector eig(it);
417 if (T(j+1,j) ==Scalar(0))
419 eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
424 eig(j) = std::complex<RealScalar>(T(j,j),T(j+1,j));
425 eig(j+1) = std::complex<RealScalar>(T(j,j+1),T(j+1,j+1));
429 if (j < it-1) eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
433 template<
typename _MatrixType,
typename _Preconditioner>
434 int DGMRES<_MatrixType, _Preconditioner>::dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, Index& neig)
const
437 typename internal::conditional<NumTraits<Scalar>::IsComplex, ComplexSchur<DenseMatrix>, RealSchur<DenseMatrix> >::type schurofH;
438 bool computeU =
true;
439 DenseMatrix matrixQ(it,it);
440 matrixQ.setIdentity();
441 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
443 ComplexVector eig(it);
444 Matrix<Index,Dynamic,1>perm(it);
445 eig = this->schurValues(schurofH);
448 DenseRealVector modulEig(it);
449 for (
int j=0; j<it; ++j) modulEig(j) = std::abs(eig(j));
450 perm.setLinSpaced(it,0,it-1);
451 internal::sortWithPermutation(modulEig, perm, neig);
455 m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
459 while (nbrEig < neig)
461 if(eig(perm(it-nbrEig-1)).imag() == RealScalar(0)) nbrEig++;
465 DenseMatrix Sr(it, nbrEig);
467 for (
int j = 0; j < nbrEig; j++)
469 Sr.col(j) = schurofH.matrixU().col(perm(it-j-1));
474 X = m_V.leftCols(it) * Sr;
478 for (
int j = 0; j < nbrEig; j++)
479 for (
int k =0; k < m_r; k++)
480 X.col(j) = X.col(j) - (m_U.col(k).dot(X.col(j)))*m_U.col(k);
484 Index m = m_V.rows();
485 if (!m_isDeflAllocated)
486 dgmresInitDeflation(m);
487 DenseMatrix MX(m, nbrEig);
489 for (
int j = 0; j < nbrEig; j++)
491 tv1 = mat * X.col(j);
492 MX.col(j) = precond.solve(tv1);
496 m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
499 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
500 m_T.block(m_r, 0, nbrEig, m_r) = X.transpose() * m_MU.leftCols(m_r);
504 for (
int j = 0; j < nbrEig; j++) m_U.col(m_r+j) = X.col(j);
505 for (
int j = 0; j < nbrEig; j++) m_MU.col(m_r+j) = MX.col(j);
510 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
513 m_isDeflInitialized =
true;
516 template<
typename _MatrixType,
typename _Preconditioner>
517 template<
typename RhsType,
typename DestType>
518 int DGMRES<_MatrixType, _Preconditioner>::dgmresApplyDeflation(
const RhsType &x, DestType &y)
const
520 DenseVector x1 = m_U.leftCols(m_r).transpose() * x;
521 y = x + m_U.leftCols(m_r) * ( m_lambdaN * m_luT.solve(x1) - x1);
527 template<
typename _MatrixType,
typename _Preconditioner,
typename Rhs>
528 struct solve_retval<DGMRES<_MatrixType, _Preconditioner>, Rhs>
529 : solve_retval_base<DGMRES<_MatrixType, _Preconditioner>, Rhs>
531 typedef DGMRES<_MatrixType, _Preconditioner> Dec;
532 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
534 template<typename Dest>
void evalTo(Dest& dst)
const
536 dec()._solve(rhs(),dst);
void setEigenv(const int neig)
Definition: DGMRES.h:196
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition: DGMRES.h:19
DGMRES(const EigenBase< MatrixDerived > &A)
Definition: DGMRES.h:137
int restart()
Definition: DGMRES.h:186
int dgmresCycle(const MatrixType &mat, const Preconditioner &precond, Dest &x, DenseVector &r0, RealScalar &beta, const RealScalar &normRhs, int &nbIts) const
Perform one restart cycle of DGMRES.
Definition: DGMRES.h:300
void setMaxEigenv(const int maxNeig)
Definition: DGMRES.h:210
DGMRES()
Definition: DGMRES.h:124
const internal::solve_retval_with_guess< DGMRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: DGMRES.h:148
int deflSize()
Definition: DGMRES.h:205
void set_restart(const int restart)
Definition: DGMRES.h:191
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:256