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>
101 class DGMRES :
public IterativeSolverBase<DGMRES<_MatrixType,_Preconditioner> >
137 template<
typename MatrixDerived>
143 template<
typename Rhs,
typename Dest>
147 for(
Index j=0; j<
b.cols(); ++j)
162 template<
typename Rhs,
typename Dest>
199 template<
typename Rhs,
typename Dest>
202 template<
typename Dest>
207 template<
typename RhsType,
typename DestType>
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);
468 Index m = m_V.rows();
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);