Warning, file /include/eigen3/unsupported/Eigen/src/IterativeSolvers/GMRES.h was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef EIGEN_GMRES_H
0012 #define EIGEN_GMRES_H
0013
0014 namespace Eigen {
0015
0016 namespace internal {
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
0056 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
0057 Index &iters, const Index &restart, typename Dest::RealScalar & tol_error) {
0058
0059 using std::sqrt;
0060 using std::abs;
0061
0062 typedef typename Dest::RealScalar RealScalar;
0063 typedef typename Dest::Scalar Scalar;
0064 typedef Matrix < Scalar, Dynamic, 1 > VectorType;
0065 typedef Matrix < Scalar, Dynamic, Dynamic, ColMajor> FMatrixType;
0066
0067 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
0068
0069 if(rhs.norm() <= considerAsZero)
0070 {
0071 x.setZero();
0072 tol_error = 0;
0073 return true;
0074 }
0075
0076 RealScalar tol = tol_error;
0077 const Index maxIters = iters;
0078 iters = 0;
0079
0080 const Index m = mat.rows();
0081
0082
0083 VectorType p0 = rhs - mat*x;
0084 VectorType r0 = precond.solve(p0);
0085
0086 const RealScalar r0Norm = r0.norm();
0087
0088
0089 if(r0Norm == 0)
0090 {
0091 tol_error = 0;
0092 return true;
0093 }
0094
0095
0096 FMatrixType H = FMatrixType::Zero(m, restart + 1);
0097 VectorType w = VectorType::Zero(restart + 1);
0098 VectorType tau = VectorType::Zero(restart + 1);
0099
0100
0101 std::vector < JacobiRotation < Scalar > > G(restart);
0102
0103
0104 VectorType t(m), v(m), workspace(m), x_new(m);
0105
0106
0107 Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
0108 RealScalar beta;
0109 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
0110 w(0) = Scalar(beta);
0111
0112 for (Index k = 1; k <= restart; ++k)
0113 {
0114 ++iters;
0115
0116 v = VectorType::Unit(m, k - 1);
0117
0118
0119
0120 for (Index i = k - 1; i >= 0; --i) {
0121 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
0122 }
0123
0124
0125 t.noalias() = mat * v;
0126 v = precond.solve(t);
0127
0128
0129
0130 for (Index i = 0; i < k; ++i) {
0131 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
0132 }
0133
0134 if (v.tail(m - k).norm() != 0.0)
0135 {
0136 if (k <= restart)
0137 {
0138
0139 Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
0140 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
0141
0142
0143 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
0144 }
0145 }
0146
0147 if (k > 1)
0148 {
0149 for (Index i = 0; i < k - 1; ++i)
0150 {
0151
0152 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
0153 }
0154 }
0155
0156 if (k<m && v(k) != (Scalar) 0)
0157 {
0158
0159 G[k - 1].makeGivens(v(k - 1), v(k));
0160
0161
0162 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
0163 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
0164 }
0165
0166
0167 H.col(k-1).head(k) = v.head(k);
0168
0169 tol_error = abs(w(k)) / r0Norm;
0170 bool stop = (k==m || tol_error < tol || iters == maxIters);
0171
0172 if (stop || k == restart)
0173 {
0174
0175 Ref<VectorType> y = w.head(k);
0176 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
0177
0178
0179 x_new.setZero();
0180 for (Index i = k - 1; i >= 0; --i)
0181 {
0182 x_new(i) += y(i);
0183
0184 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
0185 }
0186
0187 x += x_new;
0188
0189 if(stop)
0190 {
0191 return true;
0192 }
0193 else
0194 {
0195 k=0;
0196
0197
0198 p0.noalias() = rhs - mat*x;
0199 r0 = precond.solve(p0);
0200
0201
0202 H.setZero();
0203 w.setZero();
0204 tau.setZero();
0205
0206
0207 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
0208 w(0) = Scalar(beta);
0209 }
0210 }
0211 }
0212
0213 return false;
0214
0215 }
0216
0217 }
0218
0219 template< typename _MatrixType,
0220 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
0221 class GMRES;
0222
0223 namespace internal {
0224
0225 template< typename _MatrixType, typename _Preconditioner>
0226 struct traits<GMRES<_MatrixType,_Preconditioner> >
0227 {
0228 typedef _MatrixType MatrixType;
0229 typedef _Preconditioner Preconditioner;
0230 };
0231
0232 }
0233
0234
0235
0236
0237
0238
0239
0240
0241
0242
0243
0244
0245
0246
0247
0248
0249
0250
0251
0252
0253
0254
0255
0256
0257
0258
0259
0260
0261
0262
0263
0264
0265
0266
0267
0268 template< typename _MatrixType, typename _Preconditioner>
0269 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
0270 {
0271 typedef IterativeSolverBase<GMRES> Base;
0272 using Base::matrix;
0273 using Base::m_error;
0274 using Base::m_iterations;
0275 using Base::m_info;
0276 using Base::m_isInitialized;
0277
0278 private:
0279 Index m_restart;
0280
0281 public:
0282 using Base::_solve_impl;
0283 typedef _MatrixType MatrixType;
0284 typedef typename MatrixType::Scalar Scalar;
0285 typedef typename MatrixType::RealScalar RealScalar;
0286 typedef _Preconditioner Preconditioner;
0287
0288 public:
0289
0290
0291 GMRES() : Base(), m_restart(30) {}
0292
0293
0294
0295
0296
0297
0298
0299
0300
0301
0302
0303 template<typename MatrixDerived>
0304 explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
0305
0306 ~GMRES() {}
0307
0308
0309
0310 Index get_restart() { return m_restart; }
0311
0312
0313
0314
0315 void set_restart(const Index restart) { m_restart=restart; }
0316
0317
0318 template<typename Rhs,typename Dest>
0319 void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
0320 {
0321 m_iterations = Base::maxIterations();
0322 m_error = Base::m_tolerance;
0323 bool ret = internal::gmres(matrix(), b, x, Base::m_preconditioner, m_iterations, m_restart, m_error);
0324 m_info = (!ret) ? NumericalIssue
0325 : m_error <= Base::m_tolerance ? Success
0326 : NoConvergence;
0327 }
0328
0329 protected:
0330
0331 };
0332
0333 }
0334
0335 #endif