$darkmode
Eigen-unsupported  5.0.1-dev
DGMRES.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_DGMRES_H
11 #define EIGEN_DGMRES_H
12 
13 #include "../../../../Eigen/Eigenvalues"
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 template <typename MatrixType_, typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
21 class DGMRES;
22 
23 namespace internal {
24 
25 template <typename MatrixType_, typename Preconditioner_>
26 struct traits<DGMRES<MatrixType_, Preconditioner_> > {
27  typedef MatrixType_ MatrixType;
28  typedef Preconditioner_ Preconditioner;
29 };
30 
39 template <typename VectorType, typename IndexType>
40 void sortWithPermutation(VectorType& vec, IndexType& perm, typename IndexType::Scalar& ncut) {
41  eigen_assert(vec.size() == perm.size());
42  bool flag;
43  for (Index k = 0; k < ncut; k++) {
44  flag = false;
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));
48  flag = true;
49  }
50  if (!flag) break; // The vector is in sorted order
51  }
52  }
53 }
54 
55 } // namespace internal
97 template <typename MatrixType_, typename Preconditioner_>
98 class DGMRES : public IterativeSolverBase<DGMRES<MatrixType_, Preconditioner_> > {
99  typedef IterativeSolverBase<DGMRES> Base;
100  using Base::m_error;
101  using Base::m_info;
102  using Base::m_isInitialized;
103  using Base::m_iterations;
104  using Base::m_tolerance;
105  using Base::matrix;
106 
107  public:
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;
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<ComplexScalar, Dynamic, 1> ComplexVector;
121 
124  : Base(), m_restart(30), m_neig(0), m_r(0), m_maxNeig(5), m_isDeflAllocated(false), m_isDeflInitialized(false) {}
125 
136  template <typename MatrixDerived>
137  explicit DGMRES(const EigenBase<MatrixDerived>& A)
138  : Base(A.derived()),
139  m_restart(30),
140  m_neig(0),
141  m_r(0),
142  m_maxNeig(5),
143  m_isDeflAllocated(false),
144  m_isDeflInitialized(false) {}
145 
146  ~DGMRES() {}
147 
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);
153 
154  m_iterations = Base::maxIterations();
155  m_error = Base::m_tolerance;
156 
157  dgmres(matrix(), b, x, Base::m_preconditioner);
158  }
159 
163  Index restart() { return m_restart; }
164 
168  void set_restart(const Index restart) { m_restart = restart; }
169 
173  void setEigenv(const Index neig) {
174  m_neig = neig;
175  if (neig + 1 > m_maxNeig) m_maxNeig = neig + 1; // To allow for complex conjugates
176  }
177 
181  Index deflSize() { return m_r; }
182 
186  void setMaxEigenv(const Index maxNeig) { m_maxNeig = maxNeig; }
187 
188  protected:
189  // DGMRES algorithm
190  template <typename Rhs, typename Dest>
191  void dgmres(const MatrixType& mat, const Rhs& rhs, Dest& x, const Preconditioner& precond) const;
192  // Perform one cycle of GMRES
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;
196  // Compute data to use for deflation
197  Index dgmresComputeDeflationData(const MatrixType& mat, const Preconditioner& precond, const Index& it,
198  StorageIndex& neig) const;
199  // Apply deflation to a vector
200  template <typename RhsType, typename DestType>
201  Index dgmresApplyDeflation(const RhsType& In, DestType& Out) const;
202  ComplexVector schurValues(const ComplexSchur<DenseMatrix>& schurofH) const;
203  ComplexVector schurValues(const RealSchur<DenseMatrix>& schurofH) const;
204  // Init data for deflation
205  void dgmresInitDeflation(Index& rows) const;
206  mutable DenseMatrix m_V; // Krylov basis vectors
207  mutable DenseMatrix m_H; // Hessenberg matrix
208  mutable DenseMatrix m_Hes; // Initial hessenberg matrix without Givens rotations applied
209  mutable Index m_restart; // Maximum size of the Krylov subspace
210  mutable DenseMatrix m_U; // Vectors that form the basis of the invariant subspace
211  mutable DenseMatrix m_MU; // matrix operator applied to m_U (for next cycles)
212  mutable DenseMatrix m_T; /* T=U^T*M^{-1}*A*U */
213  mutable PartialPivLU<DenseMatrix> m_luT; // LU factorization of m_T
214  mutable StorageIndex m_neig; // Number of eigenvalues to extract at each restart
215  mutable Index m_r; // Current number of deflated eigenvalues, size of m_U
216  mutable Index m_maxNeig; // Maximum number of eigenvalues to deflate
217  mutable RealScalar m_lambdaN; // Modulus of the largest eigenvalue of A
218  mutable bool m_isDeflAllocated;
219  mutable bool m_isDeflInitialized;
220 
221  // Adaptive strategy
222  mutable RealScalar m_smv; // Smaller multiple of the remaining number of steps allowed
223  mutable bool m_force; // Force the use of deflation at each restart
224 };
231 template <typename MatrixType_, typename Preconditioner_>
232 template <typename Rhs, typename Dest>
233 void DGMRES<MatrixType_, Preconditioner_>::dgmres(const MatrixType& mat, const Rhs& rhs, Dest& x,
234  const Preconditioner& precond) const {
235  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
236 
237  RealScalar normRhs = rhs.norm();
238  if (normRhs <= considerAsZero) {
239  x.setZero();
240  m_error = 0;
241  return;
242  }
243 
244  // Initialization
245  m_isDeflInitialized = false;
246  Index n = mat.rows();
247  DenseVector r0(n);
248  Index nbIts = 0;
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);
252  // Initial residual vector and initial norm
253  if (x.squaredNorm() == 0) x = precond.solve(rhs);
254  r0 = rhs - mat * x;
255  RealScalar beta = r0.norm();
256 
257  m_error = beta / normRhs;
258  if (m_error < m_tolerance)
259  m_info = Success;
260  else
261  m_info = NoConvergence;
262 
263  // Iterative process
264  while (nbIts < m_iterations && m_info == NoConvergence) {
265  dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
266 
267  // Compute the new residual vector for the restart
268  if (nbIts < m_iterations && m_info == NoConvergence) {
269  r0 = rhs - mat * x;
270  beta = r0.norm();
271  }
272  }
273 }
274 
285 template <typename MatrixType_, typename Preconditioner_>
286 template <typename Dest>
287 Index DGMRES<MatrixType_, Preconditioner_>::dgmresCycle(const MatrixType& mat, const Preconditioner& precond, Dest& x,
288  DenseVector& r0, RealScalar& beta, const RealScalar& normRhs,
289  Index& nbIts) const {
290  // Initialization
291  DenseVector g(m_restart + 1); // Right hand side of the least square problem
292  g.setZero();
293  g(0) = Scalar(beta);
294  m_V.col(0) = r0 / beta;
295  m_info = NoConvergence;
296  std::vector<JacobiRotation<Scalar> > gr(m_restart); // Givens rotations
297  Index it = 0; // Number of inner iterations
298  Index n = mat.rows();
299  DenseVector tv1(n), tv2(n); // Temporary vectors
300  while (m_info == NoConvergence && it < m_restart && nbIts < m_iterations) {
301  // Apply preconditioner(s) at right
302  if (m_isDeflInitialized) {
303  dgmresApplyDeflation(m_V.col(it), tv1); // Deflation
304  tv2 = precond.solve(tv1);
305  } else {
306  tv2 = precond.solve(m_V.col(it)); // User's selected preconditioner
307  }
308  tv1 = mat * tv2;
309 
310  // Orthogonalize it with the previous basis in the basis using modified Gram-Schmidt
311  Scalar coef;
312  for (Index i = 0; i <= it; ++i) {
313  coef = tv1.dot(m_V.col(i));
314  tv1 = tv1 - coef * m_V.col(i);
315  m_H(i, it) = coef;
316  m_Hes(i, it) = coef;
317  }
318  // Normalize the vector
319  coef = tv1.norm();
320  m_V.col(it + 1) = tv1 / coef;
321  m_H(it + 1, it) = coef;
322  // m_Hes(it+1,it) = coef;
323 
324  // FIXME Check for happy breakdown
325 
326  // Update Hessenberg matrix with Givens rotations
327  for (Index i = 1; i <= it; ++i) {
328  m_H.col(it).applyOnTheLeft(i - 1, i, gr[i - 1].adjoint());
329  }
330  // Compute the new plane rotation
331  gr[it].makeGivens(m_H(it, it), m_H(it + 1, it));
332  // Apply the new rotation
333  m_H.col(it).applyOnTheLeft(it, it + 1, gr[it].adjoint());
334  g.applyOnTheLeft(it, it + 1, gr[it].adjoint());
335 
336  beta = std::abs(g(it + 1));
337  m_error = beta / normRhs;
338  // std::cerr << nbIts << " Relative Residual Norm " << m_error << std::endl;
339  it++;
340  nbIts++;
341 
342  if (m_error < m_tolerance) {
343  // The method has converged
344  m_info = Success;
345  break;
346  }
347  }
348 
349  // Compute the new coefficients by solving the least square problem
350  // it++;
351  // FIXME Check first if the matrix is singular ... zero diagonal
352  DenseVector nrs(m_restart);
353  nrs = m_H.topLeftCorner(it, it).template triangularView<Upper>().solve(g.head(it));
354 
355  // Form the new solution
356  if (m_isDeflInitialized) {
357  tv1 = m_V.leftCols(it) * nrs;
358  dgmresApplyDeflation(tv1, tv2);
359  x = x + precond.solve(tv2);
360  } else
361  x = x + precond.solve(m_V.leftCols(it) * nrs);
362 
363  // Go for a new cycle and compute data for deflation
364  if (nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r + m_neig) < m_maxNeig)
365  dgmresComputeDeflationData(mat, precond, it, m_neig);
366  return 0;
367 }
368 
369 template <typename MatrixType_, typename Preconditioner_>
371  m_U.resize(rows, m_maxNeig);
372  m_MU.resize(rows, m_maxNeig);
373  m_T.resize(m_maxNeig, m_maxNeig);
374  m_lambdaN = 0.0;
375  m_isDeflAllocated = true;
376 }
377 
378 template <typename MatrixType_, typename Preconditioner_>
379 inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
380  const ComplexSchur<DenseMatrix>& schurofH) const {
381  return schurofH.matrixT().diagonal();
382 }
383 
384 template <typename MatrixType_, typename Preconditioner_>
385 inline typename DGMRES<MatrixType_, Preconditioner_>::ComplexVector DGMRES<MatrixType_, Preconditioner_>::schurValues(
386  const RealSchur<DenseMatrix>& schurofH) const {
387  const DenseMatrix& T = schurofH.matrixT();
388  Index it = T.rows();
389  ComplexVector eig(it);
390  Index j = 0;
391  while (j < it - 1) {
392  if (T(j + 1, j) == Scalar(0)) {
393  eig(j) = ComplexScalar(T(j, j), RealScalar(0));
394  j++;
395  } else {
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));
398  j++;
399  }
400  }
401  if (j < it - 1) eig(j) = ComplexScalar(T(j, j), RealScalar(0));
402  return eig;
403 }
404 
405 template <typename MatrixType_, typename Preconditioner_>
406 Index DGMRES<MatrixType_, Preconditioner_>::dgmresComputeDeflationData(const MatrixType& mat,
407  const Preconditioner& precond, const Index& it,
408  StorageIndex& neig) const {
409  // First, find the Schur form of the Hessenberg matrix H
410  std::conditional_t<NumTraits<Scalar>::IsComplex, ComplexSchur<DenseMatrix>, RealSchur<DenseMatrix> > schurofH;
411  bool computeU = true;
412  DenseMatrix matrixQ(it, it);
413  matrixQ.setIdentity();
414  schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it, it), matrixQ, computeU);
415 
416  ComplexVector eig(it);
418  eig = this->schurValues(schurofH);
419 
420  // Reorder the absolute values of Schur values
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);
425 
426  if (!m_lambdaN) {
427  m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
428  }
429  // Count the real number of extracted eigenvalues (with complex conjugates)
430  Index nbrEig = 0;
431  while (nbrEig < neig) {
432  if (eig(perm(it - nbrEig - 1)).imag() == RealScalar(0))
433  nbrEig++;
434  else
435  nbrEig += 2;
436  }
437  // Extract the Schur vectors corresponding to the smallest Ritz values
438  DenseMatrix Sr(it, nbrEig);
439  Sr.setZero();
440  for (Index j = 0; j < nbrEig; j++) {
441  Sr.col(j) = schurofH.matrixU().col(perm(it - j - 1));
442  }
443 
444  // Form the Schur vectors of the initial matrix using the Krylov basis
445  DenseMatrix X;
446  X = m_V.leftCols(it) * Sr;
447  if (m_r) {
448  // Orthogonalize X against m_U using modified Gram-Schmidt
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);
451  }
452 
453  // Compute m_MX = A * M^-1 * X
454  Index m = m_V.rows();
455  if (!m_isDeflAllocated) dgmresInitDeflation(m);
456  DenseMatrix MX(m, nbrEig);
457  DenseVector tv1(m);
458  for (Index j = 0; j < nbrEig; j++) {
459  tv1 = mat * X.col(j);
460  MX.col(j) = precond.solve(tv1);
461  }
462 
463  // Update m_T = [U'MU U'MX; X'MU X'MX]
464  m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
465  if (m_r) {
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);
468  }
469 
470  // Save X into m_U and m_MX in m_MU
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);
473  // Increase the size of the invariant subspace
474  m_r += nbrEig;
475 
476  // Factorize m_T into m_luT
477  m_luT.compute(m_T.topLeftCorner(m_r, m_r));
478 
479  // FIXME CHeck if the factorization was correctly done (nonsingular matrix)
480  m_isDeflInitialized = true;
481  return 0;
482 }
483 template <typename MatrixType_, typename Preconditioner_>
484 template <typename RhsType, typename DestType>
485 Index 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);
488  return 0;
489 }
490 
491 } // end namespace Eigen
492 #endif
void set_restart(const Index restart)
Definition: DGMRES.h:168
DGMRES()
Definition: DGMRES.h:123
Index restart()
Definition: DGMRES.h:163
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition: DGMRES.h:21
DGMRES(const EigenBase< MatrixDerived > &A)
Definition: DGMRES.h:137
Namespace containing all symbols from the Eigen library.
void setEigenv(const Index neig)
Definition: DGMRES.h:173
Derived & setZero(Index size)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
void setMaxEigenv(const Index maxNeig)
Definition: DGMRES.h:186
Index deflSize()
Definition: DGMRES.h:181
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
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