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
18 /** \internal Low-level bi conjugate gradient stabilized algorithm
19 * \param mat The matrix A
20 * \param rhs The right hand side vector b
21 * \param x On input and initial solution, on output the computed solution.
22 * \param precond A preconditioner being able to efficiently solve for an
23 * approximation of Ax=b (regardless of b)
24 * \param iters On input the max number of iteration, on output the number of performed iterations.
25 * \param tol_error On input the tolerance error, on output an estimation of the relative error.
26 * \return false in the case of numerical issue, for example a break down of BiCGSTAB.
27 */
28 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
bicgstab(const MatrixType & mat,const Rhs & rhs,Dest & x,const Preconditioner & precond,int & iters,typename Dest::RealScalar & tol_error)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 RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
65 int i = 0;
66 int restarts = 0;
67
68 while ( r.squaredNorm()/rhs_sqnorm > tol2 && i<maxIters )
69 {
70 Scalar rho_old = rho;
71
72 rho = r0.dot(r);
73 if (abs(rho) < eps2*r0_sqnorm)
74 {
75 // The new residual vector became too orthogonal to the arbitrarily choosen direction r0
76 // Let's restart with a new r0:
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
126 /** \ingroup IterativeLinearSolvers_Module
127 * \brief A bi conjugate gradient stabilized solver for sparse square problems
128 *
129 * This class allows to solve for A.x = b sparse linear problems using a bi conjugate gradient
130 * stabilized algorithm. The vectors x and b can be either dense or sparse.
131 *
132 * \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
133 * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
134 *
135 * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
136 * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
137 * and NumTraits<Scalar>::epsilon() for the tolerance.
138 *
139 * This class can be used as the direct solver classes. Here is a typical usage example:
140 * \code
141 * int n = 10000;
142 * VectorXd x(n), b(n);
143 * SparseMatrix<double> A(n,n);
144 * // fill A and b
145 * BiCGSTAB<SparseMatrix<double> > solver;
146 * solver(A);
147 * x = solver.solve(b);
148 * std::cout << "#iterations: " << solver.iterations() << std::endl;
149 * std::cout << "estimated error: " << solver.error() << std::endl;
150 * // update b, and solve again
151 * x = solver.solve(b);
152 * \endcode
153 *
154 * By default the iterations start with x=0 as an initial guess of the solution.
155 * One can control the start using the solveWithGuess() method. Here is a step by
156 * step execution example starting with a random guess and printing the evolution
157 * of the estimated error:
158 * * \code
159 * x = VectorXd::Random(n);
160 * solver.setMaxIterations(1);
161 * int i = 0;
162 * do {
163 * x = solver.solveWithGuess(b,x);
164 * std::cout << i << " : " << solver.error() << std::endl;
165 * ++i;
166 * } while (solver.info()!=Success && i<100);
167 * \endcode
168 * Note that such a step by step excution is slightly slower.
169 *
170 * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
171 */
172 template< typename _MatrixType, typename _Preconditioner>
173 class BiCGSTAB : public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
174 {
175 typedef IterativeSolverBase<BiCGSTAB> Base;
176 using Base::mp_matrix;
177 using Base::m_error;
178 using Base::m_iterations;
179 using Base::m_info;
180 using Base::m_isInitialized;
181 public:
182 typedef _MatrixType MatrixType;
183 typedef typename MatrixType::Scalar Scalar;
184 typedef typename MatrixType::Index Index;
185 typedef typename MatrixType::RealScalar RealScalar;
186 typedef _Preconditioner Preconditioner;
187
188 public:
189
190 /** Default constructor. */
191 BiCGSTAB() : Base() {}
192
193 /** Initialize the solver with matrix \a A for further \c Ax=b solving.
194 *
195 * This constructor is a shortcut for the default constructor followed
196 * by a call to compute().
197 *
198 * \warning this class stores a reference to the matrix A as well as some
199 * precomputed values that depend on it. Therefore, if \a A is changed
200 * this class becomes invalid. Call compute() to update it with the new
201 * matrix A, or modify a copy of A.
202 */
203 BiCGSTAB(const MatrixType& A) : Base(A) {}
204
205 ~BiCGSTAB() {}
206
207 /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A
208 * \a x0 as an initial solution.
209 *
210 * \sa compute()
211 */
212 template<typename Rhs,typename Guess>
213 inline const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess>
214 solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
215 {
216 eigen_assert(m_isInitialized && "BiCGSTAB is not initialized.");
217 eigen_assert(Base::rows()==b.rows()
218 && "BiCGSTAB::solve(): invalid number of rows of the right hand side matrix b");
219 return internal::solve_retval_with_guess
220 <BiCGSTAB, Rhs, Guess>(*this, b.derived(), x0);
221 }
222
223 /** \internal */
224 template<typename Rhs,typename Dest>
225 void _solveWithGuess(const Rhs& b, Dest& x) const
226 {
227 bool failed = false;
228 for(int j=0; j<b.cols(); ++j)
229 {
230 m_iterations = Base::maxIterations();
231 m_error = Base::m_tolerance;
232
233 typename Dest::ColXpr xj(x,j);
234 if(!internal::bicgstab(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
235 failed = true;
236 }
237 m_info = failed ? NumericalIssue
238 : m_error <= Base::m_tolerance ? Success
239 : NoConvergence;
240 m_isInitialized = true;
241 }
242
243 /** \internal */
244 template<typename Rhs,typename Dest>
245 void _solve(const Rhs& b, Dest& x) const
246 {
247 // x.setZero();
248 x = b;
249 _solveWithGuess(b,x);
250 }
251
252 protected:
253
254 };
255
256
257 namespace internal {
258
259 template<typename _MatrixType, typename _Preconditioner, typename Rhs>
260 struct solve_retval<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
261 : solve_retval_base<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
262 {
263 typedef BiCGSTAB<_MatrixType, _Preconditioner> Dec;
264 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
265
266 template<typename Dest> void evalTo(Dest& dst) const
267 {
268 dec()._solve(rhs(),dst);
269 }
270 };
271
272 } // end namespace internal
273
274 } // end namespace Eigen
275
276 #endif // EIGEN_BICGSTAB_H
277