Go to the documentation of this file.
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++)
100 template<
typename _MatrixType,
typename _Preconditioner>
101 class DGMRES :
public IterativeSolverBase<DGMRES<_MatrixType,_Preconditioner> >
138 template<
typename MatrixDerived>
144 template<
typename Rhs,
typename Dest>
147 EIGEN_STATIC_ASSERT(Rhs::ColsAtCompileTime==1 || Dest::ColsAtCompileTime==1, YOU_TRIED_CALLING_A_VECTOR_METHOD_ON_A_MATRIX);
186 template<
typename Rhs,
typename Dest>
189 template<
typename Dest>
194 template<
typename RhsType,
typename DestType>
226 template<
typename _MatrixType,
typename _Preconditioner>
227 template<
typename Rhs,
typename Dest>
234 if(normRhs <= considerAsZero)
242 m_isDeflInitialized =
false;
246 m_H.resize(m_restart+1, m_restart);
247 m_Hes.resize(m_restart, m_restart);
248 m_V.resize(
n,m_restart+1);
250 if(
x.squaredNorm()==0)
251 x = precond.solve(rhs);
255 m_error =
beta/normRhs;
256 if(m_error < m_tolerance)
264 dgmresCycle(
mat, precond,
x, r0,
beta, normRhs, nbIts);
284 template<
typename _MatrixType,
typename _Preconditioner>
285 template<
typename Dest>
292 m_V.col(0) = r0/
beta;
294 std::vector<JacobiRotation<Scalar> >gr(m_restart);
298 while (m_info ==
NoConvergence && it < m_restart && nbIts < m_iterations)
301 if (m_isDeflInitialized )
303 dgmresApplyDeflation(m_V.col(it), tv1);
304 tv2 = precond.solve(tv1);
308 tv2 = precond.solve(m_V.col(it));
316 coef = tv1.dot(m_V.col(
i));
317 tv1 = tv1 - coef * m_V.col(
i);
323 m_V.col(it+1) = tv1/coef;
324 m_H(it+1, it) = coef;
332 m_H.col(it).applyOnTheLeft(
i-1,
i,gr[
i-1].
adjoint());
335 gr[it].makeGivens(m_H(it, it), m_H(it+1,it));
337 m_H.col(it).applyOnTheLeft(it,it+1,gr[it].
adjoint());
338 g.applyOnTheLeft(it,it+1, gr[it].
adjoint());
341 m_error =
beta/normRhs;
345 if (m_error < m_tolerance)
357 nrs = m_H.topLeftCorner(it,it).template triangularView<Upper>().solve(
g.head(it));
360 if (m_isDeflInitialized)
362 tv1 = m_V.leftCols(it) * nrs;
363 dgmresApplyDeflation(tv1, tv2);
364 x =
x + precond.solve(tv2);
367 x =
x + precond.solve(m_V.leftCols(it) * nrs);
370 if(nbIts < m_iterations && m_info == NoConvergence && m_neig > 0 && (m_r+m_neig) < m_maxNeig)
371 dgmresComputeDeflationData(
mat, precond, it, m_neig);
377 template<
typename _MatrixType,
typename _Preconditioner>
380 m_U.resize(
rows, m_maxNeig);
381 m_MU.resize(
rows, m_maxNeig);
382 m_T.resize(m_maxNeig, m_maxNeig);
384 m_isDeflAllocated =
true;
387 template<
typename _MatrixType,
typename _Preconditioner>
390 return schurofH.
matrixT().diagonal();
393 template<
typename _MatrixType,
typename _Preconditioner>
409 eig(
j) = std::complex<RealScalar>(
T(
j,
j),
T(
j+1,
j));
410 eig(
j+1) = std::complex<RealScalar>(
T(
j,
j+1),
T(
j+1,
j+1));
418 template<
typename _MatrixType,
typename _Preconditioner>
423 bool computeU =
true;
425 matrixQ.setIdentity();
426 schurofH.computeFromHessenberg(m_Hes.topLeftCorner(it,it), matrixQ, computeU);
430 eig = this->schurValues(schurofH);
435 perm.setLinSpaced(it,0,internal::convert_index<StorageIndex>(it-1));
440 m_lambdaN = (
std::max)(modulEig.maxCoeff(), m_lambdaN);
444 while (nbrEig < neig)
454 Sr.col(
j) = schurofH.matrixU().col(
perm(it-
j-1));
459 X = m_V.leftCols(it) * Sr;
464 for (
Index k =0; k < m_r; k++)
465 X.col(
j) =
X.col(
j) - (m_U.col(k).dot(
X.col(
j)))*m_U.col(k);
470 if (!m_isDeflAllocated)
471 dgmresInitDeflation(
m);
476 tv1 =
mat *
X.col(
j);
477 MX.col(
j) = precond.solve(tv1);
481 m_T.block(m_r, m_r, nbrEig, nbrEig) =
X.transpose() * MX;
484 m_T.block(0, m_r, m_r, nbrEig) = m_U.leftCols(m_r).transpose() * MX;
485 m_T.block(m_r, 0, nbrEig, m_r) =
X.transpose() * m_MU.leftCols(m_r);
489 for (
Index j = 0;
j < nbrEig;
j++) m_U.col(m_r+
j) =
X.col(
j);
490 for (
Index j = 0;
j < nbrEig;
j++) m_MU.col(m_r+
j) = MX.col(
j);
495 m_luT.compute(m_T.topLeftCorner(m_r, m_r));
498 m_isDeflInitialized =
true;
501 template<
typename _MatrixType,
typename _Preconditioner>
502 template<
typename RhsType,
typename DestType>
506 y =
x + m_U.leftCols(m_r) * ( m_lambdaN * m_luT.solve(
x1) -
x1);
Matrix< Scalar, Dynamic, 1 > DenseVector
void _solve_with_guess_impl(const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
idx_t idx_t idx_t idx_t idx_t * perm
Index dgmresApplyDeflation(const RhsType &In, DestType &Out) const
Namespace containing all symbols from the Eigen library.
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
LU decomposition of a matrix with partial pivoting, and related features.
void _solve_impl(const Rhs &b, Dest &x) const
IterativeSolverBase< DGMRES > Base
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Matrix< std::complex< RealScalar >, Dynamic, 1 > ComplexVector
Matrix< RealScalar, Dynamic, 1 > DenseRealVector
Eigen::Triplet< double > T
void set_restart(const Index restart)
double beta(double a, double b)
Index dgmresComputeDeflationData(const MatrixType &mat, const Preconditioner &precond, const Index &it, StorageIndex &neig) const
Performs a real Schur decomposition of a square matrix.
void dgmresInitDeflation(Index &rows) const
DGMRES(const EigenBase< MatrixDerived > &A)
MatrixType::RealScalar RealScalar
_Preconditioner Preconditioner
Index maxIterations() const
void swap(GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &a, GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &b)
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
MatrixType::Scalar Scalar
SelfAdjointEigenSolver< PlainMatrixType > eig(mat, computeVectors?ComputeEigenvectors:EigenvaluesOnly)
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.
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
Matrix< Scalar, Dynamic, 1 > DenseVector
void g(const string &key, int i)
void setMaxEigenv(const Index maxNeig)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Base class for linear iterative solvers.
NumTraits< Scalar >::Real RealScalar
MatrixType::StorageIndex StorageIndex
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Matrix< RealScalar, Dynamic, Dynamic > DenseRealMatrix
const ActualMatrixType & matrix() const
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
ComplexVector schurValues(const ComplexSchur< DenseMatrix > &schurofH) const
Preconditioner m_preconditioner
Performs a complex Schur decomposition of a real or complex square matrix.
void adjoint(const MatrixType &m)
void setEigenv(const Index neig)
PartialPivLU< DenseMatrix > m_luT
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:02:12