13#include "../../../../Eigen/Eigenvalues"
17template<
typename _MatrixType,
18 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
23template<
typename _MatrixType,
typename _Preconditioner>
38template <
typename VectorType,
typename IndexType>
46 for (
Index j = 0;
j < vec.size()-1;
j++)
48 if ( vec(perm(
j)) < vec(perm(
j+1)) )
50 std::swap(perm(
j),perm(
j+1));
100template<
typename _MatrixType,
typename _Preconditioner>
114 typedef typename MatrixType::Scalar
Scalar;
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>
226template<
typename _MatrixType,
typename _Preconditioner>
227template<
typename Rhs,
typename Dest>
231 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
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);
284template<
typename _MatrixType,
typename _Preconditioner>
285template<
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());
340 beta = std::abs(g(it+1));
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);
377template<
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;
387template<
typename _MatrixType,
typename _Preconditioner>
390 return schurofH.
matrixT().diagonal();
393template<
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));
418template<
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);
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;
501template<
typename _MatrixType,
typename _Preconditioner>
502template<
typename RhsType,
typename DestType>
506 y =
x + m_U.leftCols(m_r) * ( m_lambdaN * m_luT.solve(x1) - x1);
Matrix3f m
Definition AngleAxis_mimic_euler.cpp:1
int n
Definition BiCGSTAB_simple.cpp:1
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define eigen_assert(x)
Definition Macros.h:1037
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition StaticAssert.h:127
MatrixXf mat
Definition Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
int rows
Definition Tutorial_commainit_02.cpp:1
Eigen::Triplet< double > T
Definition Tutorial_sparse_example.cpp:6
void adjoint(const MatrixType &m)
Definition adjoint.cpp:67
Scalar * b
Definition benchVecAdd.cpp:17
SCALAR Scalar
Definition bench_gemm.cpp:46
NumTraits< Scalar >::Real RealScalar
Definition bench_gemm.cpp:47
MatrixXf MatrixType
Definition benchmark-blocking-sizes.cpp:52
Performs a complex Schur decomposition of a real or complex square matrix.
Definition ComplexSchur.h:52
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition ComplexSchur.h:162
A Restarted GMRES with deflation. This class implements a modification of the GMRES solver for sparse...
Definition DGMRES.h:102
Index m_maxNeig
Definition DGMRES.h:210
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition DGMRES.h:145
DenseMatrix m_V
Definition DGMRES.h:200
bool m_force
Definition DGMRES.h:217
DGMRES()
Definition DGMRES.h:126
DenseMatrix m_MU
Definition DGMRES.h:205
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:228
Matrix< Scalar, Dynamic, 1 > DenseVector
Definition DGMRES.h:120
RealScalar m_lambdaN
Definition DGMRES.h:211
StorageIndex m_neig
Definition DGMRES.h:208
bool m_isDeflAllocated
Definition DGMRES.h:212
DenseMatrix m_Hes
Definition DGMRES.h:202
~DGMRES()
Definition DGMRES.h:141
DenseMatrix m_U
Definition DGMRES.h:204
Matrix< RealScalar, Dynamic, Dynamic > DenseRealMatrix
Definition DGMRES.h:119
void dgmresInitDeflation(Index &rows) const
Definition DGMRES.h:378
PartialPivLU< DenseMatrix > m_luT
Definition DGMRES.h:207
Index m_r
Definition DGMRES.h:209
MatrixType::StorageIndex StorageIndex
Definition DGMRES.h:115
void setMaxEigenv(const Index maxNeig)
Definition DGMRES.h:182
Index dgmresComputeDeflationData(const MatrixType &mat, const Preconditioner &precond, const Index &it, StorageIndex &neig) const
Definition DGMRES.h:419
MatrixType::RealScalar RealScalar
Definition DGMRES.h:116
DGMRES(const EigenBase< MatrixDerived > &A)
Definition DGMRES.h:139
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition DGMRES.h:118
ComplexVector schurValues(const ComplexSchur< DenseMatrix > &schurofH) const
Definition DGMRES.h:388
Matrix< RealScalar, Dynamic, 1 > DenseRealVector
Definition DGMRES.h:121
MatrixType::Scalar Scalar
Definition DGMRES.h:114
Index m_restart
Definition DGMRES.h:203
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:286
Index deflSize()
Definition DGMRES.h:177
DenseMatrix m_H
Definition DGMRES.h:201
RealScalar m_smv
Definition DGMRES.h:216
void set_restart(const Index restart)
Definition DGMRES.h:163
_MatrixType MatrixType
Definition DGMRES.h:113
_Preconditioner Preconditioner
Definition DGMRES.h:117
DenseMatrix m_T
Definition DGMRES.h:206
void setEigenv(const Index neig)
Definition DGMRES.h:168
Matrix< std::complex< RealScalar >, Dynamic, 1 > ComplexVector
Definition DGMRES.h:122
Index dgmresApplyDeflation(const RhsType &In, DestType &Out) const
Definition DGMRES.h:503
Index restart()
Definition DGMRES.h:158
bool m_isDeflInitialized
Definition DGMRES.h:213
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:144
Index maxIterations() const
Definition IterativeSolverBase.h:281
ComputationInfo m_info
Definition IterativeSolverBase.h:438
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition IterativeSolverBase.h:250
RealScalar m_error
Definition IterativeSolverBase.h:436
void _solve_impl(const Rhs &b, Dest &x) const
Definition IterativeSolverBase.h:400
Preconditioner m_preconditioner
Definition IterativeSolverBase.h:431
void _solve_with_guess_impl(const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
Definition IterativeSolverBase.h:334
Index m_iterations
Definition IterativeSolverBase.h:437
bool m_isInitialized
Definition SparseSolverBase.h:119
DGMRES< _MatrixType, _Preconditioner > & derived()
Definition SparseSolverBase.h:79
RealScalar m_tolerance
Definition IterativeSolverBase.h:434
const ActualMatrixType & matrix() const
Definition IterativeSolverBase.h:419
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
LU decomposition of a matrix with partial pivoting, and related features.
Definition PartialPivLU.h:78
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition CwiseNullaryOp.h:562
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition PlainObjectBase.h:143
Performs a real Schur decomposition of a square matrix.
Definition RealSchur.h:55
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition RealSchur.h:144
SelfAdjointEigenSolver< PlainMatrixType > eig(mat, computeVectors?ComputeEigenvectors:EigenvaluesOnly)
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
Definition gnuplot_common_settings.hh:12
@ Success
Definition Constants.h:442
@ NoConvergence
Definition Constants.h:446
#define X
Definition icosphere.cpp:20
Scalar * y
Definition level1_cplx_impl.h:124
void sortWithPermutation(VectorType &vec, IndexType &perm, typename IndexType::Scalar &ncut)
Computes a permutation vector to have a sorted sequence.
Definition DGMRES.h:39
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
const int Dynamic
Definition Constants.h:22
Definition BandTriangularSolver.h:13
Definition EigenBase.h:30
_Preconditioner Preconditioner
Definition DGMRES.h:27
_MatrixType MatrixType
Definition DGMRES.h:26
Definition ForwardDeclarations.h:17
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2