10 #ifndef EIGEN_DGMRES_H
11 #define EIGEN_DGMRES_H
13 #include "../../../../Eigen/Eigenvalues"
19 template<
typename MatrixType_,
20 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
25 template<
typename MatrixType_,
typename Preconditioner_>
26 struct traits<DGMRES<MatrixType_,Preconditioner_> >
29 typedef Preconditioner_ Preconditioner;
40 template <
typename VectorType,
typename IndexType>
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));
102 template<
typename MatrixType_,
typename Preconditioner_>
116 typedef typename MatrixType::Scalar
Scalar;
140 template<
typename MatrixDerived>
146 template<
typename Rhs,
typename Dest>
149 EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
188 template<
typename Rhs,
typename Dest>
191 template<
typename Dest>
196 template<
typename RhsType,
typename DestType>
228 template<
typename MatrixType_,
typename Preconditioner_>
229 template<
typename Rhs,
typename Dest>
236 if(normRhs <= considerAsZero)
244 m_isDeflInitialized =
false;
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 if(
x.squaredNorm()==0)
253 x = precond.solve(rhs);
257 m_error = beta/normRhs;
258 if(m_error < m_tolerance)
266 dgmresCycle(
mat, precond,
x, r0, beta, normRhs, nbIts);
286 template<
typename MatrixType_,
typename Preconditioner_>
287 template<
typename Dest>
294 m_V.col(0) = r0/beta;
296 std::vector<JacobiRotation<Scalar> >gr(m_restart);
300 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
303 if (m_isDeflInitialized )
305 dgmresApplyDeflation(m_V.col(it), tv1);
306 tv2 = precond.solve(tv1);
310 tv2 = precond.solve(m_V.col(it));
318 coef = tv1.dot(m_V.col(
i));
319 tv1 = tv1 - coef * m_V.col(
i);
325 m_V.col(it+1) = tv1/coef;
326 m_H(it+1, it) = coef;
334 m_H.col(it).applyOnTheLeft(
i-1,
i,gr[
i-1].adjoint());
337 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
339 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].adjoint());
340 g.applyOnTheLeft(it,it+1, gr[it].adjoint());
343 m_error = beta/normRhs;
347 if (m_error < m_tolerance)
359 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(g.head(it));
362 if (m_isDeflInitialized)
364 tv1 = m_V.leftCols(it) * nrs;
365 dgmresApplyDeflation(tv1, tv2);
366 x =
x + precond.solve(tv2);
369 x =
x + precond.solve(m_V.leftCols(it) * nrs);
372 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
373 dgmresComputeDeflationData(
mat, precond, it, m_neig);
379 template<
typename MatrixType_,
typename Preconditioner_>
382 m_U.resize(
rows, m_maxNeig);
383 m_MU.resize(
rows, m_maxNeig);
384 m_T.resize(m_maxNeig, m_maxNeig);
386 m_isDeflAllocated =
true;
389 template<
typename MatrixType_,
typename Preconditioner_>
392 return schurofH.
matrixT().diagonal();
395 template<
typename MatrixType_,
typename Preconditioner_>
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));
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);
437 perm.setLinSpaced(it,0,internal::convert_index<StorageIndex>(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++;
456 Sr.col(
j) = schurofH.matrixU().col(perm(it-
j-1));
461 X = m_V.leftCols(it) * Sr;
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);
472 if (!m_isDeflAllocated)
473 dgmresInitDeflation(
m);
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 (
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);
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>
508 y =
x + m_U.leftCols(m_r) * ( m_lambdaN * m_luT.solve(x1) - x1);
SparseMatrix< double > A(n, n)
#define EIGEN_STATIC_ASSERT(X, MSG)
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...
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.
Preconditioner_ Preconditioner
void dgmresInitDeflation(Index &rows) const
DGMRES(const EigenBase< MatrixDerived > &A)
const ActualMatrixType & matrix() const
IterativeSolverBase< DGMRES > Base
Index dgmresApplyDeflation(const RhsType &In, DestType &Out) const
Matrix< Scalar, Dynamic, 1 > DenseVector
void setMaxEigenv(const Index maxNeig)
Matrix< std::complex< RealScalar >, Dynamic, 1 > ComplexVector
PartialPivLU< DenseMatrix > m_luT
Matrix< RealScalar, Dynamic, 1 > DenseRealVector
Matrix< RealScalar, Dynamic, Dynamic > DenseRealMatrix
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
void dgmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond) const
Perform several cycles of restarted GMRES with modified Gram Schmidt,.
MatrixType::StorageIndex StorageIndex
void set_restart(const Index restart)
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
MatrixType::RealScalar RealScalar
void setEigenv(const Index neig)
MatrixType::Scalar Scalar
Index dgmresComputeDeflationData(const MatrixType &mat, const Preconditioner &precond, const Index &it, StorageIndex &neig) const
ComplexVector schurValues(const ComplexSchur< DenseMatrix > &schurofH) const
TransposeReturnType transpose()
Index maxIterations() const
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
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)
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
const MatrixType & matrixT() const
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
: 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)
adouble abs(const adouble &x)