Eigen  3.2.91
BiCGSTAB.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011-2014 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, Index& 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  Index maxIters = iters;
40 
41  Index n = mat.cols();
42  VectorType r = rhs - mat * x;
43  VectorType r0 = r;
44 
45  RealScalar r0_sqnorm = r0.squaredNorm();
46  RealScalar rhs_sqnorm = rhs.squaredNorm();
47  if(rhs_sqnorm == 0)
48  {
49  x.setZero();
50  return true;
51  }
52  Scalar rho = 1;
53  Scalar alpha = 1;
54  Scalar w = 1;
55 
56  VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
57  VectorType y(n), z(n);
58  VectorType kt(n), ks(n);
59 
60  VectorType s(n), t(n);
61 
62  RealScalar tol2 = tol*tol*rhs_sqnorm;
63  RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
64  Index i = 0;
65  Index restarts = 0;
66 
67  while ( r.squaredNorm() > tol2 && i<maxIters )
68  {
69  Scalar rho_old = rho;
70 
71  rho = r0.dot(r);
72  if (abs(rho) < eps2*r0_sqnorm)
73  {
74  // The new residual vector became too orthogonal to the arbitrarily chosen direction r0
75  // Let's restart with a new r0:
76  r = rhs - mat * x;
77  r0 = r;
78  rho = r0_sqnorm = r.squaredNorm();
79  if(restarts++ == 0)
80  i = 0;
81  }
82  Scalar beta = (rho/rho_old) * (alpha / w);
83  p = r + beta * (p - w * v);
84 
85  y = precond.solve(p);
86 
87  v.noalias() = mat * y;
88 
89  alpha = rho / r0.dot(v);
90  s = r - alpha * v;
91 
92  z = precond.solve(s);
93  t.noalias() = mat * z;
94 
95  RealScalar tmp = t.squaredNorm();
96  if(tmp>RealScalar(0))
97  w = t.dot(s) / tmp;
98  else
99  w = Scalar(0);
100  x += alpha * y + w * z;
101  r = s - w * t;
102  ++i;
103  }
104  tol_error = sqrt(r.squaredNorm()/rhs_sqnorm);
105  iters = i;
106  return true;
107 }
108 
109 }
110 
111 template< typename _MatrixType,
112  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
113 class BiCGSTAB;
114 
115 namespace internal {
116 
117 template< typename _MatrixType, typename _Preconditioner>
118 struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
119 {
120  typedef _MatrixType MatrixType;
121  typedef _Preconditioner Preconditioner;
122 };
123 
124 }
125 
153 template< typename _MatrixType, typename _Preconditioner>
154 class BiCGSTAB : public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
155 {
156  typedef IterativeSolverBase<BiCGSTAB> Base;
157  using Base::mp_matrix;
158  using Base::m_error;
159  using Base::m_iterations;
160  using Base::m_info;
161  using Base::m_isInitialized;
162 public:
163  typedef _MatrixType MatrixType;
164  typedef typename MatrixType::Scalar Scalar;
165  typedef typename MatrixType::RealScalar RealScalar;
166  typedef _Preconditioner Preconditioner;
167 
168 public:
169 
171  BiCGSTAB() : Base() {}
172 
183  explicit BiCGSTAB(const MatrixType& A) : Base(A) {}
184 
185  ~BiCGSTAB() {}
186 
188  template<typename Rhs,typename Dest>
189  void _solve_with_guess_impl(const Rhs& b, Dest& x) const
190  {
191  bool failed = false;
192  for(Index j=0; j<b.cols(); ++j)
193  {
194  m_iterations = Base::maxIterations();
195  m_error = Base::m_tolerance;
196 
197  typename Dest::ColXpr xj(x,j);
198  if(!internal::bicgstab(mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
199  failed = true;
200  }
201  m_info = failed ? NumericalIssue
202  : m_error <= Base::m_tolerance ? Success
203  : NoConvergence;
204  m_isInitialized = true;
205  }
206 
208  using Base::_solve_impl;
209  template<typename Rhs,typename Dest>
210  void _solve_impl(const MatrixBase<Rhs>& b, Dest& x) const
211  {
212  x.resize(this->rows(),b.cols());
213  x.setZero();
214  _solve_with_guess_impl(b,x);
215  }
216 
217 protected:
218 
219 };
220 
221 } // end namespace Eigen
222 
223 #endif // EIGEN_BICGSTAB_H
BiCGSTAB(const MatrixType &A)
Definition: BiCGSTAB.h:183
Definition: LDLT.h:16
Index maxIterations() const
Definition: IterativeSolverBase.h:155
Definition: Constants.h:426
Definition: Constants.h:424
Definition: Eigen_Colamd.h:54
BiCGSTAB()
Definition: BiCGSTAB.h:171
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:113
Definition: Constants.h:428