DGMRES.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 template< typename MatrixType_,
20  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 {
28  typedef MatrixType_ MatrixType;
29  typedef Preconditioner_ Preconditioner;
30 };
31 
40 template <typename VectorType, typename IndexType>
41 void sortWithPermutation (VectorType& vec, IndexType& perm, typename IndexType::Scalar& ncut)
42 {
43  eigen_assert(vec.size() == perm.size());
44  bool flag;
45  for (Index k = 0; k < ncut; k++)
46  {
47  flag = false;
48  for (Index j = 0; j < vec.size()-1; j++)
49  {
50  if ( vec(perm(j)) < vec(perm(j+1)) )
51  {
52  std::swap(perm(j),perm(j+1));
53  flag = true;
54  }
55  if (!flag) break; // The vector is in sorted order
56  }
57  }
58 }
59 
60 }
102 template< typename MatrixType_, typename Preconditioner_>
103 class DGMRES : public IterativeSolverBase<DGMRES<MatrixType_,Preconditioner_> >
104 {
106  using Base::matrix;
107  using Base::m_error;
108  using Base::m_iterations;
109  using Base::m_info;
110  using Base::m_isInitialized;
111  using Base::m_tolerance;
112  public:
113  using Base::_solve_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;
125 
126 
129 
140  template<typename MatrixDerived>
142 
143  ~DGMRES() {}
144 
146  template<typename Rhs,typename Dest>
147  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
148  {
149  EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
150 
153 
155  }
156 
160  Index restart() { return m_restart; }
161 
166 
170  void setEigenv(const Index neig)
171  {
172  m_neig = neig;
173  if (neig+1 > m_maxNeig) m_maxNeig = neig+1; // To allow for complex conjugates
174  }
175 
179  Index deflSize() {return m_r; }
180 
184  void setMaxEigenv(const Index maxNeig) { m_maxNeig = maxNeig; }
185 
186  protected:
187  // DGMRES algorithm
188  template<typename Rhs, typename Dest>
189  void dgmres(const MatrixType& mat,const Rhs& rhs, Dest& x, const Preconditioner& precond) const;
190  // Perform one cycle of GMRES
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;
193  // Compute data to use for deflation
194  Index dgmresComputeDeflationData(const MatrixType& mat, const Preconditioner& precond, const Index& it, StorageIndex& neig) const;
195  // Apply deflation to a vector
196  template<typename RhsType, typename DestType>
197  Index dgmresApplyDeflation(const RhsType& In, DestType& Out) const;
198  ComplexVector schurValues(const ComplexSchur<DenseMatrix>& schurofH) const;
199  ComplexVector schurValues(const RealSchur<DenseMatrix>& schurofH) const;
200  // Init data for deflation
201  void dgmresInitDeflation(Index& rows) const;
202  mutable DenseMatrix m_V; // Krylov basis vectors
203  mutable DenseMatrix m_H; // Hessenberg matrix
204  mutable DenseMatrix m_Hes; // Initial hessenberg matrix without Givens rotations applied
205  mutable Index m_restart; // Maximum size of the Krylov subspace
206  mutable DenseMatrix m_U; // Vectors that form the basis of the invariant subspace
207  mutable DenseMatrix m_MU; // matrix operator applied to m_U (for next cycles)
208  mutable DenseMatrix m_T; /* T=U^T*M^{-1}*A*U */
209  mutable PartialPivLU<DenseMatrix> m_luT; // LU factorization of m_T
210  mutable StorageIndex m_neig; //Number of eigenvalues to extract at each restart
211  mutable Index m_r; // Current number of deflated eigenvalues, size of m_U
212  mutable Index m_maxNeig; // Maximum number of eigenvalues to deflate
213  mutable RealScalar m_lambdaN; //Modulus of the largest eigenvalue of A
214  mutable bool m_isDeflAllocated;
215  mutable bool m_isDeflInitialized;
216 
217  //Adaptive strategy
218  mutable RealScalar m_smv; // Smaller multiple of the remaining number of steps allowed
219  mutable bool m_force; // Force the use of deflation at each restart
220 
221 };
228 template< typename MatrixType_, typename Preconditioner_>
229 template<typename Rhs, typename Dest>
231  const Preconditioner& precond) const
232 {
233  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
234 
235  RealScalar normRhs = rhs.norm();
236  if(normRhs <= considerAsZero)
237  {
238  x.setZero();
239  m_error = 0;
240  return;
241  }
242 
243  //Initialization
244  m_isDeflInitialized = false;
245  Index n = mat.rows();
246  DenseVector r0(n);
247  Index nbIts = 0;
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);
251  //Initial residual vector and initial norm
252  if(x.squaredNorm()==0)
253  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  {
266  dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
267 
268  // Compute the new residual vector for the restart
269  if (nbIts < m_iterations && m_info == NoConvergence) {
270  r0 = rhs - mat * x;
271  beta = r0.norm();
272  }
273  }
274 }
275 
286 template< typename MatrixType_, typename Preconditioner_>
287 template<typename Dest>
288 Index DGMRES<MatrixType_, Preconditioner_>::dgmresCycle(const MatrixType& mat, const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta, const RealScalar& normRhs, Index& nbIts) const
289 {
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  {
302  // Apply preconditioner(s) at right
303  if (m_isDeflInitialized )
304  {
305  dgmresApplyDeflation(m_V.col(it), tv1); // Deflation
306  tv2 = precond.solve(tv1);
307  }
308  else
309  {
310  tv2 = precond.solve(m_V.col(it)); // User's selected preconditioner
311  }
312  tv1 = mat * tv2;
313 
314  // Orthogonalize it with the previous basis in the basis using modified Gram-Schmidt
315  Scalar coef;
316  for (Index i = 0; i <= it; ++i)
317  {
318  coef = tv1.dot(m_V.col(i));
319  tv1 = tv1 - coef * m_V.col(i);
320  m_H(i,it) = coef;
321  m_Hes(i,it) = coef;
322  }
323  // Normalize the vector
324  coef = tv1.norm();
325  m_V.col(it+1) = tv1/coef;
326  m_H(it+1, it) = coef;
327 // m_Hes(it+1,it) = coef;
328 
329  // FIXME Check for happy breakdown
330 
331  // Update Hessenberg matrix with Givens rotations
332  for (Index i = 1; i <= it; ++i)
333  {
334  m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
335  }
336  // Compute the new plane rotation
337  gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
338  // Apply the new rotation
339  m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
340  g.applyOnTheLeft(it,it+1, gr[it].adjoint());
341 
342  beta = std::abs(g(it+1));
343  m_error = beta/normRhs;
344  // std::cerr << nbIts << " Relative Residual Norm " << m_error << std::endl;
345  it++; nbIts++;
346 
347  if (m_error < m_tolerance)
348  {
349  // The method has converged
350  m_info = Success;
351  break;
352  }
353  }
354 
355  // Compute the new coefficients by solving the least square problem
356 // it++;
357  //FIXME Check first if the matrix is singular ... zero diagonal
358  DenseVector nrs(m_restart);
359  nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
360 
361  // Form the new solution
362  if (m_isDeflInitialized)
363  {
364  tv1 = m_V.leftCols(it) * nrs;
365  dgmresApplyDeflation(tv1, tv2);
366  x = x + precond.solve(tv2);
367  }
368  else
369  x = x + precond.solve(m_V.leftCols(it) * nrs);
370 
371  // Go for a new cycle and compute data for deflation
372  if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
373  dgmresComputeDeflationData(mat, precond, it, m_neig);
374  return 0;
375 
376 }
377 
378 
379 template< typename MatrixType_, typename Preconditioner_>
381 {
382  m_U.resize(rows, m_maxNeig);
383  m_MU.resize(rows, m_maxNeig);
384  m_T.resize(m_maxNeig, m_maxNeig);
385  m_lambdaN = 0.0;
386  m_isDeflAllocated = true;
387 }
388 
389 template< typename MatrixType_, typename Preconditioner_>
391 {
392  return schurofH.matrixT().diagonal();
393 }
394 
395 template< typename MatrixType_, typename Preconditioner_>
397 {
398  const DenseMatrix& T = schurofH.matrixT();
399  Index it = T.rows();
400  ComplexVector eig(it);
401  Index j = 0;
402  while (j < it-1)
403  {
404  if (T(j+1,j) ==Scalar(0))
405  {
406  eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
407  j++;
408  }
409  else
410  {
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));
413  j++;
414  }
415  }
416  if (j < it-1) eig(j) = std::complex<RealScalar>(T(j,j),RealScalar(0));
417  return eig;
418 }
419 
420 template< typename MatrixType_, typename Preconditioner_>
422 {
423  // First, find the Schur form of the Hessenberg matrix H
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);
429 
430  ComplexVector eig(it);
432  eig = this->schurValues(schurofH);
433 
434  // Reorder the absolute values of Schur values
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);
439 
440  if (!m_lambdaN)
441  {
442  m_lambdaN = (std::max)(modulEig.maxCoeff(), m_lambdaN);
443  }
444  //Count the real number of extracted eigenvalues (with complex conjugates)
445  Index nbrEig = 0;
446  while (nbrEig < neig)
447  {
448  if(eig(perm(it-nbrEig-1)).imag() == RealScalar(0)) nbrEig++;
449  else nbrEig += 2;
450  }
451  // Extract the Schur vectors corresponding to the smallest Ritz values
452  DenseMatrix Sr(it, nbrEig);
453  Sr.setZero();
454  for (Index j = 0; j < nbrEig; j++)
455  {
456  Sr.col(j) = schurofH.matrixU().col(perm(it-j-1));
457  }
458 
459  // Form the Schur vectors of the initial matrix using the Krylov basis
460  DenseMatrix X;
461  X = m_V.leftCols(it) * Sr;
462  if (m_r)
463  {
464  // Orthogonalize X against m_U using modified Gram-Schmidt
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);
468  }
469 
470  // Compute m_MX = A * M^-1 * X
471  Index m = m_V.rows();
472  if (!m_isDeflAllocated)
473  dgmresInitDeflation(m);
474  DenseMatrix MX(m, nbrEig);
475  DenseVector tv1(m);
476  for (Index j = 0; j < nbrEig; j++)
477  {
478  tv1 = mat * X.col(j);
479  MX.col(j) = precond.solve(tv1);
480  }
481 
482  //Update m_T = [U'MU U'MX; X'MU X'MX]
483  m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
484  if(m_r)
485  {
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);
488  }
489 
490  // Save X into m_U and m_MX in m_MU
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);
493  // Increase the size of the invariant subspace
494  m_r += nbrEig;
495 
496  // Factorize m_T into m_luT
497  m_luT.compute(m_T.topLeftCorner(m_r, m_r));
498 
499  //FIXME CHeck if the factorization was correctly done (nonsingular matrix)
500  m_isDeflInitialized = true;
501  return 0;
502 }
503 template<typename MatrixType_, typename Preconditioner_>
504 template<typename RhsType, typename DestType>
506 {
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);
509  return 0;
510 }
511 
512 } // end namespace Eigen
513 #endif
Matrix3f m
int n
SparseMatrix< double > A(n, n)
int i
Matrix3f y
#define eigen_assert(x)
MatrixXf X
#define EIGEN_STATIC_ASSERT(X, MSG)
MatrixXf mat
Matrix< float, 1, Dynamic > MatrixType
const ComplexMatrixType & matrixT() const
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition: DGMRES.h:104
RealScalar m_error
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
Preconditioner_ Preconditioner
Definition: DGMRES.h:119
Index restart()
Definition: DGMRES.h:160
void dgmresInitDeflation(Index &rows) const
Definition: DGMRES.h:380
DGMRES(const EigenBase< MatrixDerived > &A)
Definition: DGMRES.h:141
const ActualMatrixType & matrix() const
RealScalar m_smv
Definition: DGMRES.h:218
IterativeSolverBase< DGMRES > Base
Definition: DGMRES.h:105
Index dgmresApplyDeflation(const RhsType &In, DestType &Out) const
Definition: DGMRES.h:505
Matrix< Scalar, Dynamic, 1 > DenseVector
Definition: DGMRES.h:122
void setMaxEigenv(const Index maxNeig)
Definition: DGMRES.h:184
DenseMatrix m_H
Definition: DGMRES.h:203
Matrix< std::complex< RealScalar >, Dynamic, 1 > ComplexVector
Definition: DGMRES.h:124
DenseMatrix m_V
Definition: DGMRES.h:202
PartialPivLU< DenseMatrix > m_luT
Definition: DGMRES.h:209
Matrix< RealScalar, Dynamic, 1 > DenseRealVector
Definition: DGMRES.h:123
Matrix< RealScalar, Dynamic, Dynamic > DenseRealMatrix
Definition: DGMRES.h:121
StorageIndex m_neig
Definition: DGMRES.h:210
DenseMatrix m_T
Definition: DGMRES.h:208
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: DGMRES.h:147
Index m_iterations
DenseMatrix m_U
Definition: DGMRES.h:206
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
bool m_isDeflAllocated
Definition: DGMRES.h:214
MatrixType::StorageIndex StorageIndex
Definition: DGMRES.h:117
void set_restart(const Index restart)
Definition: DGMRES.h:165
Index m_restart
Definition: DGMRES.h:205
DenseMatrix m_MU
Definition: DGMRES.h:207
DenseMatrix m_Hes
Definition: DGMRES.h:204
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: DGMRES.h:120
bool m_force
Definition: DGMRES.h:219
MatrixType::RealScalar RealScalar
Definition: DGMRES.h:118
Index deflSize()
Definition: DGMRES.h:179
void setEigenv(const Index neig)
Definition: DGMRES.h:170
MatrixType::Scalar Scalar
Definition: DGMRES.h:116
Index m_maxNeig
Definition: DGMRES.h:212
Index dgmresComputeDeflationData(const MatrixType &mat, const Preconditioner &precond, const Index &it, StorageIndex &neig) const
Definition: DGMRES.h:421
Index m_r
Definition: DGMRES.h:211
bool m_isDeflInitialized
Definition: DGMRES.h:215
ComplexVector schurValues(const ComplexSchur< DenseMatrix > &schurofH) const
Definition: DGMRES.h:390
MatrixType_ MatrixType
Definition: DGMRES.h:115
RealScalar m_lambdaN
Definition: DGMRES.h:213
TransposeReturnType transpose()
std::enable_if_t< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 > _solve_with_guess_impl(const Rhs &b, MatrixBase< DestDerived > &aDest) const
void _solve_impl(const Rhs &b, Dest &x) const
Preconditioner m_preconditioner
const ActualMatrixType & matrix() const
Derived & setZero(Index rows, Index cols)
const MatrixType & matrixT() const
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
Definition: DGMRES.h:41
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
CleanedUpDerType< DerType >::type() max(const AutoDiffScalar< DerType > &x, const T &y)
const int Dynamic
adouble abs(const adouble &x)
Definition: AdolcForward:74
Derived::Index rows
std::ptrdiff_t j