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


shape_reconstruction
Author(s): Roberto Martín-Martín
autogenerated on Sat Jun 8 2019 18:31:29