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>
38 template <
typename VectorType,
typename IndexType>
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)) )
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;
125 DGMRES() : Base(),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
137 template<
typename MatrixDerived>
138 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) {}
143 template<
typename Rhs,
typename Dest>
147 for(
Index j=0; j<b.cols(); ++j)
149 m_iterations = Base::maxIterations();
150 m_error = Base::m_tolerance;
153 dgmres(
matrix(), b.col(j), xj, Base::m_preconditioner);
156 : m_error <= Base::m_tolerance ?
Success 158 m_isInitialized =
true;
162 template<
typename Rhs,
typename Dest>
166 _solve_with_guess_impl(b,x.derived());
184 if (neig+1 > m_maxNeig) m_maxNeig = neig+1;
199 template<
typename Rhs,
typename Dest>
200 void dgmres(
const MatrixType&
mat,
const Rhs& rhs, Dest&
x,
const Preconditioner& precond)
const;
202 template<
typename Dest>
203 Index dgmresCycle(
const MatrixType& mat,
const Preconditioner& precond, Dest& x, DenseVector& r0, RealScalar& beta,
const RealScalar& normRhs,
Index& nbIts)
const;
205 Index dgmresComputeDeflationData(
const MatrixType& mat,
const Preconditioner& precond,
const Index& it, StorageIndex& neig)
const;
207 template<
typename RhsType,
typename DestType>
208 Index dgmresApplyDeflation(
const RhsType& In, DestType& Out)
const;
212 void dgmresInitDeflation(
Index& rows)
const;
239 template<
typename _MatrixType,
typename _Preconditioner>
240 template<
typename Rhs,
typename Dest>
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);
252 x = precond.solve(x);
256 m_error = beta/normRhs;
257 if(m_error < m_tolerance)
265 dgmresCycle(mat, precond, x, r0, beta, normRhs, nbIts);
283 template<
typename _MatrixType,
typename _Preconditioner>
284 template<
typename Dest>
291 m_V.col(0) = r0/beta;
293 std::vector<JacobiRotation<Scalar> >gr(m_restart);
297 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
300 if (m_isDeflInitialized )
302 dgmresApplyDeflation(m_V.col(it), tv1);
303 tv2 = precond.solve(tv1);
307 tv2 = precond.solve(m_V.col(it));
313 for (
Index i = 0; i <= it; ++i)
315 coef = tv1.dot(m_V.col(i));
316 tv1 = tv1 - coef * m_V.col(i);
322 m_V.col(it+1) = tv1/coef;
323 m_H(it+1, it) = coef;
329 for (
Index i = 1; i <= it; ++i)
331 m_H.col(it).applyOnTheLeft(i-1,i,gr[i-1].adjoint());
334 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
336 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
337 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
340 m_error = beta/normRhs;
344 if (m_error < m_tolerance)
356 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
359 if (m_isDeflInitialized)
361 tv1 = m_V.leftCols(it) * nrs;
362 dgmresApplyDeflation(tv1, tv2);
363 x = x + precond.solve(tv2);
366 x = x + precond.solve(m_V.leftCols(it) * nrs);
369 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
370 dgmresComputeDeflationData(mat, precond, it, m_neig);
376 template<
typename _MatrixType,
typename _Preconditioner>
379 m_U.resize(rows, m_maxNeig);
380 m_MU.resize(rows, m_maxNeig);
381 m_T.resize(m_maxNeig, m_maxNeig);
383 m_isDeflAllocated =
true;
386 template<
typename _MatrixType,
typename _Preconditioner>
389 return schurofH.
matrixT().diagonal();
392 template<
typename _MatrixType,
typename _Preconditioner>
408 eig(j) = std::complex<RealScalar>(T(j,j),T(j+1,j));
409 eig(j+1) = std::complex<RealScalar>(T(j,j+1),T(j+1,j+1));
413 if (j < it-1)
eig(j) = std::complex<RealScalar>(T(j,j),
RealScalar(0));
417 template<
typename _MatrixType,
typename _Preconditioner>
422 bool computeU =
true;
424 matrixQ.setIdentity();
425 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
429 eig = this->schurValues(schurofH);
434 perm.setLinSpaced(it,0,internal::convert_index<StorageIndex>(it-1));
439 m_lambdaN = (
std::max)(modulEig.maxCoeff(), m_lambdaN);
443 while (nbrEig < neig)
451 for (
Index j = 0; j < nbrEig; j++)
453 Sr.col(j) = schurofH.matrixU().col(perm(it-j-1));
458 X = m_V.leftCols(it) * Sr;
462 for (
Index j = 0; j < nbrEig; j++)
463 for (
Index k =0; k < m_r; k++)
464 X.col(j) = X.col(j) - (m_U.col(k).dot(X.col(j)))*m_U.col(k);
469 if (!m_isDeflAllocated)
470 dgmresInitDeflation(m);
473 for (
Index j = 0; j < nbrEig; j++)
475 tv1 = mat * X.col(j);
476 MX.col(j) = precond.solve(tv1);
480 m_T.block(m_r, m_r, nbrEig, nbrEig) = X.transpose() * MX;
483 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
484 m_T.block(m_r, 0, nbrEig, m_r) = X.transpose() * m_MU.leftCols(m_r);
488 for (
Index j = 0; j < nbrEig; j++) m_U.col(m_r+j) = X.col(j);
489 for (
Index j = 0; j < nbrEig; j++) m_MU.col(m_r+j) = MX.col(j);
494 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
497 m_isDeflInitialized =
true;
500 template<
typename _MatrixType,
typename _Preconditioner>
501 template<
typename RhsType,
typename DestType>
505 y = x + m_U.leftCols(m_r) * ( m_lambdaN * m_luT.solve(x1) - x1);
Index dgmresApplyDeflation(const RhsType &In, DestType &Out) const
void _solve_with_guess_impl(const Rhs &b, Dest &x) const
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
ComplexVector schurValues(const ComplexSchur< DenseMatrix > &schurofH) const
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Performs a real Schur decomposition of a square matrix.
void dgmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond) const
Perform several cycles of restarted GMRES with modified Gram Schmidt,.
_Preconditioner Preconditioner
void setEigenv(const Index neig)
LU decomposition of a matrix with partial pivoting, and related features.
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1, !IsRowMajor > ColXpr
MatrixType::RealScalar RealScalar
Index dgmresComputeDeflationData(const MatrixType &mat, const Preconditioner &precond, const Index &it, StorageIndex &neig) const
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
MatrixType A(a, *n, *n, *lda)
MatrixType::Scalar Scalar
NumTraits< Scalar >::Real RealScalar
int EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
void set_restart(const Index restart)
Matrix< RealScalar, Dynamic, 1 > DenseRealVector
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
PartialPivLU< DenseMatrix > m_luT
SelfAdjointEigenSolver< PlainMatrixType > eig(mat, computeVectors?ComputeEigenvectors:EigenvaluesOnly)
Matrix< RealScalar, Dynamic, Dynamic > DenseRealMatrix
IterativeSolverBase< DGMRES > Base
_Preconditioner Preconditioner
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.
Matrix< std::complex< RealScalar >, Dynamic, 1 > ComplexVector
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
void _solve_impl(const Rhs &b, MatrixBase< Dest > &x) const
MatrixType::StorageIndex StorageIndex
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
void setMaxEigenv(const Index maxNeig)
Performs a complex Schur decomposition of a real or complex square matrix.
Base class for linear iterative solvers.
Base class for all dense matrices, vectors, and expressions.
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
void dgmresInitDeflation(Index &rows) const
Matrix< Scalar, Dynamic, 1 > DenseVector
DGMRES(const EigenBase< MatrixDerived > &A)