BiCGSTAB.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) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@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 #ifndef EIGEN_BICGSTAB_H
12 #define EIGEN_BICGSTAB_H
13 
14 namespace Eigen {
15 
16 namespace internal {
17 
28 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
29 bool bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x,
30  const Preconditioner& precond, int& iters,
31  typename Dest::RealScalar& tol_error)
32 {
33  using std::sqrt;
34  using std::abs;
35  typedef typename Dest::RealScalar RealScalar;
36  typedef typename Dest::Scalar Scalar;
37  typedef Matrix<Scalar,Dynamic,1> VectorType;
38  RealScalar tol = tol_error;
39  int maxIters = iters;
40 
41  int n = mat.cols();
42  x = precond.solve(x);
43  VectorType r = rhs - mat * x;
44  VectorType r0 = r;
45 
46  RealScalar r0_sqnorm = r0.squaredNorm();
47  RealScalar rhs_sqnorm = rhs.squaredNorm();
48  if(rhs_sqnorm == 0)
49  {
50  x.setZero();
51  return true;
52  }
53  Scalar rho = 1;
54  Scalar alpha = 1;
55  Scalar w = 1;
56 
57  VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
58  VectorType y(n), z(n);
59  VectorType kt(n), ks(n);
60 
61  VectorType s(n), t(n);
62 
63  RealScalar tol2 = tol*tol;
64  int i = 0;
65  int restarts = 0;
66 
67  while ( r.squaredNorm()/rhs_sqnorm > tol2 && i<maxIters )
68  {
69  Scalar rho_old = rho;
70 
71  rho = r0.dot(r);
72  if (internal::isMuchSmallerThan(rho,r0_sqnorm))
73  {
74  // The new residual vector became too orthogonal to the arbitrarily choosen direction r0
75  // Let's restart with a new r0:
76  r0 = r;
77  rho = r0_sqnorm = r.squaredNorm();
78  if(restarts++ == 0)
79  i = 0;
80  }
81  Scalar beta = (rho/rho_old) * (alpha / w);
82  p = r + beta * (p - w * v);
83 
84  y = precond.solve(p);
85 
86  v.noalias() = mat * y;
87 
88  alpha = rho / r0.dot(v);
89  s = r - alpha * v;
90 
91  z = precond.solve(s);
92  t.noalias() = mat * z;
93 
94  RealScalar tmp = t.squaredNorm();
95  if(tmp>RealScalar(0))
96  w = t.dot(s) / tmp;
97  else
98  w = Scalar(0);
99  x += alpha * y + w * z;
100  r = s - w * t;
101  ++i;
102  }
103  tol_error = sqrt(r.squaredNorm()/rhs_sqnorm);
104  iters = i;
105  return true;
106 }
107 
108 }
109 
110 template< typename _MatrixType,
111  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
112 class BiCGSTAB;
113 
114 namespace internal {
115 
116 template< typename _MatrixType, typename _Preconditioner>
117 struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
118 {
119  typedef _MatrixType MatrixType;
120  typedef _Preconditioner Preconditioner;
121 };
122 
123 }
124 
171 template< typename _MatrixType, typename _Preconditioner>
172 class BiCGSTAB : public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
173 {
175  using Base::mp_matrix;
176  using Base::m_error;
177  using Base::m_iterations;
178  using Base::m_info;
179  using Base::m_isInitialized;
180 public:
181  typedef _MatrixType MatrixType;
182  typedef typename MatrixType::Scalar Scalar;
183  typedef typename MatrixType::Index Index;
184  typedef typename MatrixType::RealScalar RealScalar;
185  typedef _Preconditioner Preconditioner;
186 
187 public:
188 
190  BiCGSTAB() : Base() {}
191 
202  BiCGSTAB(const MatrixType& A) : Base(A) {}
203 
205 
211  template<typename Rhs,typename Guess>
213  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
214  {
215  eigen_assert(m_isInitialized && "BiCGSTAB is not initialized.");
216  eigen_assert(Base::rows()==b.rows()
217  && "BiCGSTAB::solve(): invalid number of rows of the right hand side matrix b");
219  <BiCGSTAB, Rhs, Guess>(*this, b.derived(), x0);
220  }
221 
223  template<typename Rhs,typename Dest>
224  void _solveWithGuess(const Rhs& b, Dest& x) const
225  {
226  bool failed = false;
227  for(int j=0; j<b.cols(); ++j)
228  {
229  m_iterations = Base::maxIterations();
230  m_error = Base::m_tolerance;
231 
232  typename Dest::ColXpr xj(x,j);
233  if(!internal::bicgstab(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
234  failed = true;
235  }
236  m_info = failed ? NumericalIssue
237  : m_error <= Base::m_tolerance ? Success
238  : NoConvergence;
239  m_isInitialized = true;
240  }
241 
243  template<typename Rhs,typename Dest>
244  void _solve(const Rhs& b, Dest& x) const
245  {
246 // x.setZero();
247  x = b;
248  _solveWithGuess(b,x);
249  }
250 
251 protected:
252 
253 };
254 
255 
256 namespace internal {
257 
258  template<typename _MatrixType, typename _Preconditioner, typename Rhs>
259 struct solve_retval<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
260  : solve_retval_base<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
261 {
263  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
264 
265  template<typename Dest> void evalTo(Dest& dst) const
266  {
267  dec()._solve(rhs(),dst);
268  }
269 };
270 
271 } // end namespace internal
272 
273 } // end namespace Eigen
274 
275 #endif // EIGEN_BICGSTAB_H
A preconditioner based on the digonal entries.
void _solve(const Rhs &b, Dest &x) const
Definition: BiCGSTAB.h:244
IntermediateState sqrt(const Expression &arg)
_MatrixType MatrixType
Definition: BiCGSTAB.h:181
BiCGSTAB(const MatrixType &A)
Definition: BiCGSTAB.h:202
_Preconditioner Preconditioner
Definition: BiCGSTAB.h:185
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
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, typename NumTraits< Scalar >::Real precision=NumTraits< Scalar >::dummy_precision())
void _solveWithGuess(const Rhs &b, Dest &x) const
Definition: BiCGSTAB.h:224
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
const internal::solve_retval_with_guess< BiCGSTAB, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: BiCGSTAB.h:213
MatrixType::Index Index
Definition: BiCGSTAB.h:183
#define v
void rhs(const real_t *x, real_t *f)
MatrixType::RealScalar RealScalar
Definition: BiCGSTAB.h:184
bool bicgstab(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, int &iters, typename Dest::RealScalar &tol_error)
Definition: BiCGSTAB.h:29
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
Definition: Solve.h:61
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:112
MatrixType::Scalar Scalar
Definition: BiCGSTAB.h:182
IterativeSolverBase< BiCGSTAB > Base
Definition: BiCGSTAB.h:174
EIGEN_STRONG_INLINE Index cols() const
#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:28