MINRES.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Giacomo Po <gpo@ucla.edu>
5 // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 
12 #ifndef EIGEN_MINRES_H_
13 #define EIGEN_MINRES_H_
14 
15 
16 namespace Eigen {
17 
18  namespace internal {
19 
29  template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
31  void minres(const MatrixType& mat, const Rhs& rhs, Dest& x,
32  const Preconditioner& precond, int& iters,
33  typename Dest::RealScalar& tol_error)
34  {
35  using std::sqrt;
36  typedef typename Dest::RealScalar RealScalar;
37  typedef typename Dest::Scalar Scalar;
38  typedef Matrix<Scalar,Dynamic,1> VectorType;
39 
40  // initialize
41  const int maxIters(iters); // initialize maxIters to iters
42  const int N(mat.cols()); // the size of the matrix
43  const RealScalar rhsNorm2(rhs.squaredNorm());
44  const RealScalar threshold2(tol_error*tol_error*rhsNorm2); // convergence threshold (compared to residualNorm2)
45 
46  // Initialize preconditioned Lanczos
47 // VectorType v_old(N); // will be initialized inside loop
48  VectorType v( VectorType::Zero(N) ); //initialize v
49  VectorType v_new(rhs-mat*x); //initialize v_new
50  RealScalar residualNorm2(v_new.squaredNorm());
51 // VectorType w(N); // will be initialized inside loop
52  VectorType w_new(precond.solve(v_new)); // initialize w_new
53 // RealScalar beta; // will be initialized inside loop
54  RealScalar beta_new2(v_new.dot(w_new));
55  eigen_assert(beta_new2 >= 0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
56  RealScalar beta_new(sqrt(beta_new2));
57  const RealScalar beta_one(beta_new);
58  v_new /= beta_new;
59  w_new /= beta_new;
60  // Initialize other variables
61  RealScalar c(1.0); // the cosine of the Givens rotation
62  RealScalar c_old(1.0);
63  RealScalar s(0.0); // the sine of the Givens rotation
64  RealScalar s_old(0.0); // the sine of the Givens rotation
65 // VectorType p_oold(N); // will be initialized in loop
66  VectorType p_old(VectorType::Zero(N)); // initialize p_old=0
67  VectorType p(p_old); // initialize p=0
68  RealScalar eta(1.0);
69 
70  iters = 0; // reset iters
71  while ( iters < maxIters ){
72 
73  // Preconditioned Lanczos
74  /* Note that there are 4 variants on the Lanczos algorithm. These are
75  * described in Paige, C. C. (1972). Computational variants of
76  * the Lanczos method for the eigenproblem. IMA Journal of Applied
77  * Mathematics, 10(3), 373–381. The current implementation corresponds
78  * to the case A(2,7) in the paper. It also corresponds to
79  * algorithm 6.14 in Y. Saad, Iterative Methods for Sparse Linear
80  * Systems, 2003 p.173. For the preconditioned version see
81  * A. Greenbaum, Iterative Methods for Solving Linear Systems, SIAM (1987).
82  */
83  const RealScalar beta(beta_new);
84 // v_old = v; // update: at first time step, this makes v_old = 0 so value of beta doesn't matter
85  const VectorType v_old(v); // NOT SURE IF CREATING v_old EVERY ITERATION IS EFFICIENT
86  v = v_new; // update
87 // w = w_new; // update
88  const VectorType w(w_new); // NOT SURE IF CREATING w EVERY ITERATION IS EFFICIENT
89  v_new.noalias() = mat*w - beta*v_old; // compute v_new
90  const RealScalar alpha = v_new.dot(w);
91  v_new -= alpha*v; // overwrite v_new
92  w_new = precond.solve(v_new); // overwrite w_new
93  beta_new2 = v_new.dot(w_new); // compute beta_new
94  eigen_assert(beta_new2 >= 0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
95  beta_new = sqrt(beta_new2); // compute beta_new
96  v_new /= beta_new; // overwrite v_new for next iteration
97  w_new /= beta_new; // overwrite w_new for next iteration
98 
99  // Givens rotation
100  const RealScalar r2 =s*alpha+c*c_old*beta; // s, s_old, c and c_old are still from previous iteration
101  const RealScalar r3 =s_old*beta; // s, s_old, c and c_old are still from previous iteration
102  const RealScalar r1_hat=c*alpha-c_old*s*beta;
103  const RealScalar r1 =sqrt( std::pow(r1_hat,2) + std::pow(beta_new,2) );
104  c_old = c; // store for next iteration
105  s_old = s; // store for next iteration
106  c=r1_hat/r1; // new cosine
107  s=beta_new/r1; // new sine
108 
109  // Update solution
110 // p_oold = p_old;
111  const VectorType p_oold(p_old); // NOT SURE IF CREATING p_oold EVERY ITERATION IS EFFICIENT
112  p_old = p;
113  p.noalias()=(w-r2*p_old-r3*p_oold) /r1; // IS NOALIAS REQUIRED?
114  x += beta_one*c*eta*p;
115  residualNorm2 *= s*s;
116 
117  if ( residualNorm2 < threshold2){
118  break;
119  }
120 
121  eta=-s*eta; // update eta
122  iters++; // increment iteration number (for output purposes)
123  }
124  tol_error = std::sqrt(residualNorm2 / rhsNorm2); // return error. Note that this is the estimated error. The real error |Ax-b|/|b| may be slightly larger
125  }
126 
127  }
128 
129  template< typename _MatrixType, int _UpLo=Lower,
130  typename _Preconditioner = IdentityPreconditioner>
131 // typename _Preconditioner = IdentityPreconditioner<typename _MatrixType::Scalar> > // preconditioner must be positive definite
132  class MINRES;
133 
134  namespace internal {
135 
136  template< typename _MatrixType, int _UpLo, typename _Preconditioner>
137  struct traits<MINRES<_MatrixType,_UpLo,_Preconditioner> >
138  {
139  typedef _MatrixType MatrixType;
140  typedef _Preconditioner Preconditioner;
141  };
142 
143  }
144 
194  template< typename _MatrixType, int _UpLo, typename _Preconditioner>
195  class MINRES : public IterativeSolverBase<MINRES<_MatrixType,_UpLo,_Preconditioner> >
196  {
197 
199  using Base::mp_matrix;
200  using Base::m_error;
201  using Base::m_iterations;
202  using Base::m_info;
203  using Base::m_isInitialized;
204  public:
205  typedef _MatrixType MatrixType;
206  typedef typename MatrixType::Scalar Scalar;
207  typedef typename MatrixType::Index Index;
208  typedef typename MatrixType::RealScalar RealScalar;
209  typedef _Preconditioner Preconditioner;
210 
211  enum {UpLo = _UpLo};
212 
213  public:
214 
216  MINRES() : Base() {}
217 
228  MINRES(const MatrixType& A) : Base(A) {}
229 
232 
238  template<typename Rhs,typename Guess>
240  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
241  {
242  eigen_assert(m_isInitialized && "MINRES is not initialized.");
243  eigen_assert(Base::rows()==b.rows()
244  && "MINRES::solve(): invalid number of rows of the right hand side matrix b");
246  <MINRES, Rhs, Guess>(*this, b.derived(), x0);
247  }
248 
250  template<typename Rhs,typename Dest>
251  void _solveWithGuess(const Rhs& b, Dest& x) const
252  {
253  m_iterations = Base::maxIterations();
254  m_error = Base::m_tolerance;
255 
256  for(int j=0; j<b.cols(); ++j)
257  {
258  m_iterations = Base::maxIterations();
259  m_error = Base::m_tolerance;
260 
261  typename Dest::ColXpr xj(x,j);
262  internal::minres(mp_matrix->template selfadjointView<UpLo>(), b.col(j), xj,
263  Base::m_preconditioner, m_iterations, m_error);
264  }
265 
266  m_isInitialized = true;
267  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
268  }
269 
271  template<typename Rhs,typename Dest>
272  void _solve(const Rhs& b, Dest& x) const
273  {
274  x.setZero();
275  _solveWithGuess(b,x);
276  }
277 
278  protected:
279 
280  };
281 
282  namespace internal {
283 
284  template<typename _MatrixType, int _UpLo, typename _Preconditioner, typename Rhs>
285  struct solve_retval<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
286  : solve_retval_base<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
287  {
289  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
290 
291  template<typename Dest> void evalTo(Dest& dst) const
292  {
293  dec()._solve(rhs(),dst);
294  }
295  };
296 
297  } // end namespace internal
298 
299 } // end namespace Eigen
300 
301 #endif // EIGEN_MINRES_H
302 
MatrixType::Scalar Scalar
Definition: MINRES.h:206
#define N
EIGEN_DONT_INLINE void minres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, int &iters, typename Dest::RealScalar &tol_error)
Definition: MINRES.h:31
IntermediateState sqrt(const Expression &arg)
void _solve(const Rhs &b, Dest &x) const
Definition: MINRES.h:272
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
Definition: BlockMethods.h:15
IntermediateState pow(const Expression &arg1, const Expression &arg2)
_MatrixType MatrixType
Definition: MINRES.h:205
MatrixType::Index Index
Definition: MINRES.h:207
const internal::solve_retval_with_guess< MINRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: MINRES.h:240
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:132
MINRES(const MatrixType &A)
Definition: MINRES.h:228
void _solveWithGuess(const Rhs &b, Dest &x) const
Definition: MINRES.h:251
IterativeSolverBase< MINRES > Base
Definition: MINRES.h:198
#define v
void rhs(const real_t *x, real_t *f)
_Preconditioner Preconditioner
Definition: MINRES.h:209
A naive preconditioner which approximates any matrix as the identity matrix.
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
Definition: Solve.h:61
MatrixType::RealScalar RealScalar
Definition: MINRES.h:208
#define EIGEN_DONT_INLINE
#define eigen_assert(x)
Base class for linear iterative solvers.
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:54