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 < Scalar, Dynamic, 1 > VectorType;
00065         typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
00066 
00067         RealScalar tol = tol_error;
00068         const int maxIters = iters;
00069         iters = 0;
00070 
00071         const int m = mat.rows();
00072 
00073         VectorType p0 = rhs - mat*x;
00074         VectorType r0 = precond.solve(p0);
00075 //      RealScalar r0_sqnorm = r0.squaredNorm();
00076 
00077         VectorType w = VectorType::Zero(restart + 1);
00078 
00079         FMatrixType H = FMatrixType::Zero(m, restart + 1);
00080         VectorType tau = VectorType::Zero(restart + 1);
00081         std::vector < JacobiRotation < Scalar > > G(restart);
00082 
00083         // generate first Householder vector
00084         VectorType e;
00085         RealScalar beta;
00086         r0.makeHouseholder(e, tau.coeffRef(0), beta);
00087         w(0)=(Scalar) beta;
00088         H.bottomLeftCorner(m - 1, 1) = e;
00089 
00090         for (int k = 1; k <= restart; ++k) {
00091 
00092                 ++iters;
00093 
00094                 VectorType v = VectorType::Unit(m, k - 1), workspace(m);
00095 
00096                 // apply Householder reflections H_{1} ... H_{k-1} to v
00097                 for (int i = k - 1; i >= 0; --i) {
00098                         v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00099                 }
00100 
00101                 // apply matrix M to v:  v = mat * v;
00102                 VectorType t=mat*v;
00103                 v=precond.solve(t);
00104 
00105                 // apply Householder reflections H_{k-1} ... H_{1} to v
00106                 for (int i = 0; i < k; ++i) {
00107                         v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00108                 }
00109 
00110                 if (v.tail(m - k).norm() != 0.0) {
00111 
00112                         if (k <= restart) {
00113 
00114                                 // generate new Householder vector
00115                                   VectorType e(m - k - 1);
00116                                 RealScalar beta;
00117                                 v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
00118                                 H.col(k).tail(m - k - 1) = e;
00119 
00120                                 // apply Householder reflection H_{k} to v
00121                                 v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
00122 
00123                         }
00124                 }
00125 
00126                 if (k > 1) {
00127                         for (int i = 0; i < k - 1; ++i) {
00128                                 // apply old Givens rotations to v
00129                                 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
00130                         }
00131                 }
00132 
00133                 if (k<m && v(k) != (Scalar) 0) {
00134                         // determine next Givens rotation
00135                         G[k - 1].makeGivens(v(k - 1), v(k));
00136 
00137                         // apply Givens rotation to v and w
00138                         v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00139                         w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
00140 
00141                 }
00142 
00143                 // insert coefficients into upper matrix triangle
00144                 H.col(k - 1).head(k) = v.head(k);
00145 
00146                 bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
00147 
00148                 if (stop || k == restart) {
00149 
00150                         // solve upper triangular system
00151                         VectorType y = w.head(k);
00152                         H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
00153 
00154                         // use Horner-like scheme to calculate solution vector
00155                         VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
00156 
00157                         // apply Householder reflection H_{k} to x_new
00158                         x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
00159 
00160                         for (int i = k - 2; i >= 0; --i) {
00161                                 x_new += y(i) * VectorType::Unit(m, i);
00162                                 // apply Householder reflection H_{i} to x_new
00163                                 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
00164                         }
00165 
00166                         x += x_new;
00167 
00168                         if (stop) {
00169                                 return true;
00170                         } else {
00171                                 k=0;
00172 
00173                                 // reset data for a restart  r0 = rhs - mat * x;
00174                                 VectorType p0=mat*x;
00175                                 VectorType p1=precond.solve(p0);
00176                                 r0 = rhs - p1;
00177 //                                 r0_sqnorm = r0.squaredNorm();
00178                                 w = VectorType::Zero(restart + 1);
00179                                 H = FMatrixType::Zero(m, restart + 1);
00180                                 tau = VectorType::Zero(restart + 1);
00181 
00182                                 // generate first Householder vector
00183                                 RealScalar beta;
00184                                 r0.makeHouseholder(e, tau.coeffRef(0), beta);
00185                                 w(0)=(Scalar) beta;
00186                                 H.bottomLeftCorner(m - 1, 1) = e;
00187 
00188                         }
00189 
00190                 }
00191 
00192 
00193 
00194         }
00195         
00196         return false;
00197 
00198 }
00199 
00200 }
00201 
00202 template< typename _MatrixType,
00203           typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
00204 class GMRES;
00205 
00206 namespace internal {
00207 
00208 template< typename _MatrixType, typename _Preconditioner>
00209 struct traits<GMRES<_MatrixType,_Preconditioner> >
00210 {
00211   typedef _MatrixType MatrixType;
00212   typedef _Preconditioner Preconditioner;
00213 };
00214 
00215 }
00216 
00262 template< typename _MatrixType, typename _Preconditioner>
00263 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
00264 {
00265   typedef IterativeSolverBase<GMRES> Base;
00266   using Base::mp_matrix;
00267   using Base::m_error;
00268   using Base::m_iterations;
00269   using Base::m_info;
00270   using Base::m_isInitialized;
00271  
00272 private:
00273   int m_restart;
00274   
00275 public:
00276   typedef _MatrixType MatrixType;
00277   typedef typename MatrixType::Scalar Scalar;
00278   typedef typename MatrixType::Index Index;
00279   typedef typename MatrixType::RealScalar RealScalar;
00280   typedef _Preconditioner Preconditioner;
00281 
00282 public:
00283 
00285   GMRES() : Base(), m_restart(30) {}
00286 
00297   GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
00298 
00299   ~GMRES() {}
00300   
00303   int get_restart() { return m_restart; }
00304   
00308   void set_restart(const int restart) { m_restart=restart; }
00309   
00315   template<typename Rhs,typename Guess>
00316   inline const internal::solve_retval_with_guess<GMRES, Rhs, Guess>
00317   solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
00318   {
00319     eigen_assert(m_isInitialized && "GMRES is not initialized.");
00320     eigen_assert(Base::rows()==b.rows()
00321               && "GMRES::solve(): invalid number of rows of the right hand side matrix b");
00322     return internal::solve_retval_with_guess
00323             <GMRES, Rhs, Guess>(*this, b.derived(), x0);
00324   }
00325   
00327   template<typename Rhs,typename Dest>
00328   void _solveWithGuess(const Rhs& b, Dest& x) const
00329   {    
00330     bool failed = false;
00331     for(int j=0; j<b.cols(); ++j)
00332     {
00333       m_iterations = Base::maxIterations();
00334       m_error = Base::m_tolerance;
00335       
00336       typename Dest::ColXpr xj(x,j);
00337       if(!internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
00338         failed = true;
00339     }
00340     m_info = failed ? NumericalIssue
00341            : m_error <= Base::m_tolerance ? Success
00342            : NoConvergence;
00343     m_isInitialized = true;
00344   }
00345 
00347   template<typename Rhs,typename Dest>
00348   void _solve(const Rhs& b, Dest& x) const
00349   {
00350     x = b;
00351     if(x.squaredNorm() == 0) return; // Check Zero right hand side
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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Sat Jun 8 2019 19:37:16