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>
26 typedef _MatrixType MatrixType;
27 typedef _Preconditioner Preconditioner;
38 template <
typename VectorType,
typename IndexType>
41 eigen_assert(vec.size() == perm.size());
42 typedef typename IndexType::Scalar
Index;
44 for (Index k = 0; k < ncut; k++)
47 for (Index j = 0; j < vec.size()-1; j++)
49 if ( vec(perm(j)) < vec(perm(j+1)) )
51 std::swap(perm(j),perm(j+1));
101 template<
typename _MatrixType,
typename _Preconditioner>
107 using Base::m_iterations;
109 using Base::m_isInitialized;
110 using Base::m_tolerance;
112 using Base::_solve_impl;
113 typedef _MatrixType MatrixType;
114 typedef typename MatrixType::Scalar Scalar;
115 typedef typename MatrixType::Index
Index;
116 typedef typename MatrixType::StorageIndex StorageIndex;
117 typedef typename MatrixType::RealScalar RealScalar;
118 typedef _Preconditioner Preconditioner;
127 DGMRES() : Base(),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
139 template<
typename MatrixDerived>
140 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) {}
145 template<
typename Rhs,
typename Dest>
146 void _solve_with_guess_impl(
const Rhs& b, Dest& x)
const 149 for(
int j=0; j<b.cols(); ++j)
151 m_iterations = Base::maxIterations();
152 m_error = Base::m_tolerance;
154 typename Dest::ColXpr xj(x,j);
155 dgmres(matrix(), b.col(j), xj, Base::m_preconditioner);
158 : m_error <= Base::m_tolerance ?
Success 160 m_isInitialized =
true;
164 template<
typename Rhs,
typename Dest>
168 _solve_with_guess_impl(b,x.derived());
186 if (neig+1 > m_maxNeig) m_maxNeig = neig+1;
201 template<
typename Rhs,
typename Dest>
202 void dgmres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond)
const;
204 template<
typename Dest>
205 int dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
int& nbIts)
const;
207 int dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, StorageIndex& neig)
const;
209 template<
typename RhsType,
typename DestType>
210 int dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
214 void dgmresInitDeflation(Index& rows)
const;
215 mutable DenseMatrix m_V;
216 mutable DenseMatrix m_H;
217 mutable DenseMatrix m_Hes;
218 mutable Index m_restart;
219 mutable DenseMatrix m_U;
220 mutable DenseMatrix m_MU;
221 mutable DenseMatrix m_T;
223 mutable StorageIndex m_neig;
225 mutable int m_maxNeig;
226 mutable RealScalar m_lambdaN;
227 mutable bool m_isDeflAllocated;
228 mutable bool m_isDeflInitialized;
231 mutable RealScalar m_smv;
232 mutable bool m_force;
241 template<
typename _MatrixType,
typename _Preconditioner>
242 template<
typename Rhs,
typename Dest>
244 const Preconditioner& precond)
const 250 m_H.resize(m_restart+1, m_restart);
251 m_Hes.resize(m_restart, m_restart);
252 m_V.resize(n,m_restart+1);
254 x = precond.solve(x);
256 RealScalar beta = r0.norm();
257 RealScalar normRhs = rhs.norm();
258 m_error = beta/normRhs;
259 if(m_error < m_tolerance)
267 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
285 template<
typename _MatrixType,
typename _Preconditioner>
286 template<
typename Dest>
293 m_V.col(0) = r0/beta;
295 std::vector<JacobiRotation<Scalar> >gr(m_restart);
299 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
302 if (m_isDeflInitialized )
304 dgmresApplyDeflation(m_V.col(it), tv1);
305 tv2 = precond.solve(tv1);
309 tv2 = precond.solve(m_V.col(it));
315 for (
int i = 0; i <= it; ++i)
317 coef = tv1.dot(m_V.col(i));
318 tv1 = tv1 - coef * m_V.col(i);
324 m_V.col(it+1) = tv1/coef;
325 m_H(it+1, it) = coef;
331 for (
int i = 1; i <= it; ++i)
333 m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
336 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
338 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
339 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
341 beta = std::abs(g(it+1));
342 m_error = beta/normRhs;
346 if (m_error < m_tolerance)
358 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
361 if (m_isDeflInitialized)
363 tv1 = m_V.leftCols(it) * nrs;
364 dgmresApplyDeflation(tv1, tv2);
365 x = x + precond.solve(tv2);
368 x = x + precond.solve(m_V.leftCols(it) * nrs);
371 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
372 dgmresComputeDeflationData(mat, precond, it, m_neig);
378 template<
typename _MatrixType,
typename _Preconditioner>
381 m_U.resize(rows, m_maxNeig);
382 m_MU.resize(rows, m_maxNeig);
383 m_T.resize(m_maxNeig, m_maxNeig);
385 m_isDeflAllocated =
true;
388 template<
typename _MatrixType,
typename _Preconditioner>
391 return schurofH.
matrixT().diagonal();
394 template<
typename _MatrixType,
typename _Preconditioner>
397 typedef typename MatrixType::Index
Index;
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>
425 bool computeU =
true;
427 matrixQ.setIdentity();
428 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
432 eig = this->schurValues(schurofH);
436 for (
int j=0; j<it; ++j) modulEig(j) = std::abs(eig(j));
437 perm.setLinSpaced(it,0,it-1);
442 m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
446 while (nbrEig < neig)
448 if(eig(perm(it-nbrEig-1)).imag() == RealScalar(0)) nbrEig++;
454 for (
int 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 (
int j = 0; j < nbrEig; j++)
466 for (
int 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);
476 for (
int 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 (
int j = 0; j < nbrEig; j++) m_U.col(m_r+j) = X.col(j);
492 for (
int 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>
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);
int restart()
Get the restart value.
Definition: DGMRES.h:173
A Restarted GMRES with deflation.
Definition: DGMRES.h:19
Definition: RealSchur.h:54
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:243
int deflSize()
Get the size of the deflation subspace size.
Definition: DGMRES.h:192
LU decomposition of a matrix with partial pivoting, and related features.
Definition: ForwardDeclarations.h:250
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
Definition: DGMRES.h:39
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
The provided data did not satisfy the prerequisites.
Definition: Constants.h:434
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Resizes to the given size, and sets all coefficients in this expression to zero.
Definition: CwiseNullaryOp.h:515
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
Computation was successful.
Definition: Constants.h:432
DGMRES()
Default constructor.
Definition: DGMRES.h:127
Definition: BandTriangularSolver.h:13
void set_restart(const int restart)
Set the restart value (default is 30)
Definition: DGMRES.h:178
void setEigenv(const int neig)
Set the number of eigenvalues to deflate at each restart.
Definition: DGMRES.h:183
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
Definition: ComplexSchur.h:51
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:143
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:287
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
void setMaxEigenv(const int maxNeig)
Set the maximum size of the deflation subspace.
Definition: DGMRES.h:197
Iterative procedure did not converge.
Definition: Constants.h:436
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition: ComplexSchur.h:162
Definition: ForwardDeclarations.h:17
DGMRES(const EigenBase< MatrixDerived > &A)
Initialize the solver with matrix A for further Ax=b solving.
Definition: DGMRES.h:140