00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifndef EIGEN_GMRES_H
00027 #define EIGEN_GMRES_H
00028
00029 namespace Eigen {
00030
00031 namespace internal {
00032
00070 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
00071 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
00072 int &iters, const int &restart, typename Dest::RealScalar & tol_error) {
00073
00074 using std::sqrt;
00075 using std::abs;
00076
00077 typedef typename Dest::RealScalar RealScalar;
00078 typedef typename Dest::Scalar Scalar;
00079 typedef Matrix < RealScalar, Dynamic, 1 > RealVectorType;
00080 typedef Matrix < Scalar, Dynamic, 1 > VectorType;
00081 typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
00082
00083 RealScalar tol = tol_error;
00084 const int maxIters = iters;
00085 iters = 0;
00086
00087 const int m = mat.rows();
00088
00089 VectorType p0 = rhs - mat*x;
00090 VectorType r0 = precond.solve(p0);
00091
00092
00093 VectorType w = VectorType::Zero(restart + 1);
00094
00095 FMatrixType H = FMatrixType::Zero(m, restart + 1);
00096 VectorType tau = VectorType::Zero(restart + 1);
00097 std::vector < JacobiRotation < Scalar > > G(restart);
00098
00099
00100 VectorType e;
00101 RealScalar beta;
00102 r0.makeHouseholder(e, tau.coeffRef(0), beta);
00103 w(0)=(Scalar) beta;
00104 H.bottomLeftCorner(m - 1, 1) = e;
00105
00106 for (int k = 1; k <= restart; ++k) {
00107
00108 ++iters;
00109
00110 VectorType v = VectorType::Unit(m, k - 1), workspace(m);
00111
00112
00113 for (int i = k - 1; i >= 0; --i) {
00114 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00115 }
00116
00117
00118 VectorType t=mat*v;
00119 v=precond.solve(t);
00120
00121
00122 for (int i = 0; i < k; ++i) {
00123 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00124 }
00125
00126 if (v.tail(m - k).norm() != 0.0) {
00127
00128 if (k <= restart) {
00129
00130
00131 VectorType e;
00132 RealScalar beta;
00133 v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
00134 H.col(k).tail(m - k - 1) = e;
00135
00136
00137 v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
00138
00139 }
00140 }
00141
00142 if (k > 1) {
00143 for (int i = 0; i < k - 1; ++i) {
00144
00145 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
00146 }
00147 }
00148
00149 if (k<m && v(k) != (Scalar) 0) {
00150
00151 G[k - 1].makeGivens(v(k - 1), v(k));
00152
00153
00154 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00155 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00156
00157 }
00158
00159
00160 H.col(k - 1).head(k) = v.head(k);
00161
00162 bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
00163
00164 if (stop || k == restart) {
00165
00166
00167 VectorType y = w.head(k);
00168 H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
00169
00170
00171 VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
00172
00173
00174 x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
00175
00176 for (int i = k - 2; i >= 0; --i) {
00177 x_new += y(i) * VectorType::Unit(m, i);
00178
00179 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00180 }
00181
00182 x += x_new;
00183
00184 if (stop) {
00185 return true;
00186 } else {
00187 k=0;
00188
00189
00190 VectorType p0=mat*x;
00191 VectorType p1=precond.solve(p0);
00192 r0 = rhs - p1;
00193
00194 w = VectorType::Zero(restart + 1);
00195 H = FMatrixType::Zero(m, restart + 1);
00196 tau = VectorType::Zero(restart + 1);
00197
00198
00199 RealScalar beta;
00200 r0.makeHouseholder(e, tau.coeffRef(0), beta);
00201 w(0)=(Scalar) beta;
00202 H.bottomLeftCorner(m - 1, 1) = e;
00203
00204 }
00205
00206 }
00207
00208
00209
00210 }
00211
00212 return false;
00213
00214 }
00215
00216 }
00217
00218 template< typename _MatrixType,
00219 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
00220 class GMRES;
00221
00222 namespace internal {
00223
00224 template< typename _MatrixType, typename _Preconditioner>
00225 struct traits<GMRES<_MatrixType,_Preconditioner> >
00226 {
00227 typedef _MatrixType MatrixType;
00228 typedef _Preconditioner Preconditioner;
00229 };
00230
00231 }
00232
00278 template< typename _MatrixType, typename _Preconditioner>
00279 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
00280 {
00281 typedef IterativeSolverBase<GMRES> Base;
00282 using Base::mp_matrix;
00283 using Base::m_error;
00284 using Base::m_iterations;
00285 using Base::m_info;
00286 using Base::m_isInitialized;
00287
00288 private:
00289 int m_restart;
00290
00291 public:
00292 typedef _MatrixType MatrixType;
00293 typedef typename MatrixType::Scalar Scalar;
00294 typedef typename MatrixType::Index Index;
00295 typedef typename MatrixType::RealScalar RealScalar;
00296 typedef _Preconditioner Preconditioner;
00297
00298 public:
00299
00301 GMRES() : Base(), m_restart(30) {}
00302
00313 GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
00314
00315 ~GMRES() {}
00316
00319 int get_restart() { return m_restart; }
00320
00324 void set_restart(const int restart) { m_restart=restart; }
00325
00331 template<typename Rhs,typename Guess>
00332 inline const internal::solve_retval_with_guess<GMRES, Rhs, Guess>
00333 solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
00334 {
00335 eigen_assert(m_isInitialized && "GMRES is not initialized.");
00336 eigen_assert(Base::rows()==b.rows()
00337 && "GMRES::solve(): invalid number of rows of the right hand side matrix b");
00338 return internal::solve_retval_with_guess
00339 <GMRES, Rhs, Guess>(*this, b.derived(), x0);
00340 }
00341
00343 template<typename Rhs,typename Dest>
00344 void _solveWithGuess(const Rhs& b, Dest& x) const
00345 {
00346 bool failed = false;
00347 for(int j=0; j<b.cols(); ++j)
00348 {
00349 m_iterations = Base::maxIterations();
00350 m_error = Base::m_tolerance;
00351
00352 typename Dest::ColXpr xj(x,j);
00353 if(!internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
00354 failed = true;
00355 }
00356 m_info = failed ? NumericalIssue
00357 : m_error <= Base::m_tolerance ? Success
00358 : NoConvergence;
00359 m_isInitialized = true;
00360 }
00361
00363 template<typename Rhs,typename Dest>
00364 void _solve(const Rhs& b, Dest& x) const
00365 {
00366 x.setZero();
00367 _solveWithGuess(b,x);
00368 }
00369
00370 protected:
00371
00372 };
00373
00374
00375 namespace internal {
00376
00377 template<typename _MatrixType, typename _Preconditioner, typename Rhs>
00378 struct solve_retval<GMRES<_MatrixType, _Preconditioner>, Rhs>
00379 : solve_retval_base<GMRES<_MatrixType, _Preconditioner>, Rhs>
00380 {
00381 typedef GMRES<_MatrixType, _Preconditioner> Dec;
00382 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
00383
00384 template<typename Dest> void evalTo(Dest& dst) const
00385 {
00386 dec()._solve(rhs(),dst);
00387 }
00388 };
00389
00390 }
00391
00392 }
00393
00394 #endif // EIGEN_GMRES_H