GMRES.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2012 Kolja Brix <brix@igpm.rwth-aaachen.de>
00006 //
00007 // This Source Code Form is subject to the terms of the Mozilla
00008 // Public License v. 2.0. If a copy of the MPL was not distributed
00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00010 
00011 #ifndef EIGEN_GMRES_H
00012 #define EIGEN_GMRES_H
00013 
00014 namespace Eigen { 
00015 
00016 namespace internal {
00017 
00055 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
00056 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
00057                 int &iters, const int &restart, typename Dest::RealScalar & tol_error) {
00058 
00059         using std::sqrt;
00060         using std::abs;
00061 
00062         typedef typename Dest::RealScalar RealScalar;
00063         typedef typename Dest::Scalar Scalar;
00064         typedef Matrix < RealScalar, Dynamic, 1 > RealVectorType;
00065         typedef Matrix < Scalar, Dynamic, 1 > VectorType;
00066         typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
00067 
00068         RealScalar tol = tol_error;
00069         const int maxIters = iters;
00070         iters = 0;
00071 
00072         const int m = mat.rows();
00073 
00074         VectorType p0 = rhs - mat*x;
00075         VectorType r0 = precond.solve(p0);
00076 //      RealScalar r0_sqnorm = r0.squaredNorm();
00077 
00078         VectorType w = VectorType::Zero(restart + 1);
00079 
00080         FMatrixType H = FMatrixType::Zero(m, restart + 1);
00081         VectorType tau = VectorType::Zero(restart + 1);
00082         std::vector < JacobiRotation < Scalar > > G(restart);
00083 
00084         // generate first Householder vector
00085         VectorType e;
00086         RealScalar beta;
00087         r0.makeHouseholder(e, tau.coeffRef(0), beta);
00088         w(0)=(Scalar) beta;
00089         H.bottomLeftCorner(m - 1, 1) = e;
00090 
00091         for (int k = 1; k <= restart; ++k) {
00092 
00093                 ++iters;
00094 
00095                 VectorType v = VectorType::Unit(m, k - 1), workspace(m);
00096 
00097                 // apply Householder reflections H_{1} ... H_{k-1} to v
00098                 for (int i = k - 1; i >= 0; --i) {
00099                         v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00100                 }
00101 
00102                 // apply matrix M to v:  v = mat * v;
00103                 VectorType t=mat*v;
00104                 v=precond.solve(t);
00105 
00106                 // apply Householder reflections H_{k-1} ... H_{1} to v
00107                 for (int i = 0; i < k; ++i) {
00108                         v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00109                 }
00110 
00111                 if (v.tail(m - k).norm() != 0.0) {
00112 
00113                         if (k <= restart) {
00114 
00115                                 // generate new Householder vector
00116                                   VectorType e(m - k - 1);
00117                                 RealScalar beta;
00118                                 v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
00119                                 H.col(k).tail(m - k - 1) = e;
00120 
00121                                 // apply Householder reflection H_{k} to v
00122                                 v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
00123 
00124                         }
00125                 }
00126 
00127                 if (k > 1) {
00128                         for (int i = 0; i < k - 1; ++i) {
00129                                 // apply old Givens rotations to v
00130                                 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
00131                         }
00132                 }
00133 
00134                 if (k<m && v(k) != (Scalar) 0) {
00135                         // determine next Givens rotation
00136                         G[k - 1].makeGivens(v(k - 1), v(k));
00137 
00138                         // apply Givens rotation to v and w
00139                         v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00140                         w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00141 
00142                 }
00143 
00144                 // insert coefficients into upper matrix triangle
00145                 H.col(k - 1).head(k) = v.head(k);
00146 
00147                 bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
00148 
00149                 if (stop || k == restart) {
00150 
00151                         // solve upper triangular system
00152                         VectorType y = w.head(k);
00153                         H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
00154 
00155                         // use Horner-like scheme to calculate solution vector
00156                         VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
00157 
00158                         // apply Householder reflection H_{k} to x_new
00159                         x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
00160 
00161                         for (int i = k - 2; i >= 0; --i) {
00162                                 x_new += y(i) * VectorType::Unit(m, i);
00163                                 // apply Householder reflection H_{i} to x_new
00164                                 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00165                         }
00166 
00167                         x += x_new;
00168 
00169                         if (stop) {
00170                                 return true;
00171                         } else {
00172                                 k=0;
00173 
00174                                 // reset data for a restart  r0 = rhs - mat * x;
00175                                 VectorType p0=mat*x;
00176                                 VectorType p1=precond.solve(p0);
00177                                 r0 = rhs - p1;
00178 //                                 r0_sqnorm = r0.squaredNorm();
00179                                 w = VectorType::Zero(restart + 1);
00180                                 H = FMatrixType::Zero(m, restart + 1);
00181                                 tau = VectorType::Zero(restart + 1);
00182 
00183                                 // generate first Householder vector
00184                                 RealScalar beta;
00185                                 r0.makeHouseholder(e, tau.coeffRef(0), beta);
00186                                 w(0)=(Scalar) beta;
00187                                 H.bottomLeftCorner(m - 1, 1) = e;
00188 
00189                         }
00190 
00191                 }
00192 
00193 
00194 
00195         }
00196         
00197         return false;
00198 
00199 }
00200 
00201 }
00202 
00203 template< typename _MatrixType,
00204           typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
00205 class GMRES;
00206 
00207 namespace internal {
00208 
00209 template< typename _MatrixType, typename _Preconditioner>
00210 struct traits<GMRES<_MatrixType,_Preconditioner> >
00211 {
00212   typedef _MatrixType MatrixType;
00213   typedef _Preconditioner Preconditioner;
00214 };
00215 
00216 }
00217 
00263 template< typename _MatrixType, typename _Preconditioner>
00264 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
00265 {
00266   typedef IterativeSolverBase<GMRES> Base;
00267   using Base::mp_matrix;
00268   using Base::m_error;
00269   using Base::m_iterations;
00270   using Base::m_info;
00271   using Base::m_isInitialized;
00272  
00273 private:
00274   int m_restart;
00275   
00276 public:
00277   typedef _MatrixType MatrixType;
00278   typedef typename MatrixType::Scalar Scalar;
00279   typedef typename MatrixType::Index Index;
00280   typedef typename MatrixType::RealScalar RealScalar;
00281   typedef _Preconditioner Preconditioner;
00282 
00283 public:
00284 
00286   GMRES() : Base(), m_restart(30) {}
00287 
00298   GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
00299 
00300   ~GMRES() {}
00301   
00304   int get_restart() { return m_restart; }
00305   
00309   void set_restart(const int restart) { m_restart=restart; }
00310   
00316   template<typename Rhs,typename Guess>
00317   inline const internal::solve_retval_with_guess<GMRES, Rhs, Guess>
00318   solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
00319   {
00320     eigen_assert(m_isInitialized && "GMRES is not initialized.");
00321     eigen_assert(Base::rows()==b.rows()
00322               && "GMRES::solve(): invalid number of rows of the right hand side matrix b");
00323     return internal::solve_retval_with_guess
00324             <GMRES, Rhs, Guess>(*this, b.derived(), x0);
00325   }
00326   
00328   template<typename Rhs,typename Dest>
00329   void _solveWithGuess(const Rhs& b, Dest& x) const
00330   {    
00331     bool failed = false;
00332     for(int j=0; j<b.cols(); ++j)
00333     {
00334       m_iterations = Base::maxIterations();
00335       m_error = Base::m_tolerance;
00336       
00337       typename Dest::ColXpr xj(x,j);
00338       if(!internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
00339         failed = true;
00340     }
00341     m_info = failed ? NumericalIssue
00342            : m_error <= Base::m_tolerance ? Success
00343            : NoConvergence;
00344     m_isInitialized = true;
00345   }
00346 
00348   template<typename Rhs,typename Dest>
00349   void _solve(const Rhs& b, Dest& x) const
00350   {
00351     x.setZero();
00352     _solveWithGuess(b,x);
00353   }
00354 
00355 protected:
00356 
00357 };
00358 
00359 
00360 namespace internal {
00361 
00362   template<typename _MatrixType, typename _Preconditioner, typename Rhs>
00363 struct solve_retval<GMRES<_MatrixType, _Preconditioner>, Rhs>
00364   : solve_retval_base<GMRES<_MatrixType, _Preconditioner>, Rhs>
00365 {
00366   typedef GMRES<_MatrixType, _Preconditioner> Dec;
00367   EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
00368 
00369   template<typename Dest> void evalTo(Dest& dst) const
00370   {
00371     dec()._solve(rhs(),dst);
00372   }
00373 };
00374 
00375 } // end namespace internal
00376 
00377 } // end namespace Eigen
00378 
00379 #endif // EIGEN_GMRES_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Mon Oct 6 2014 12:24:35