11 #ifndef EIGEN_BICGSTAB_H
12 #define EIGEN_BICGSTAB_H
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)
35 typedef typename Dest::RealScalar RealScalar;
36 typedef typename Dest::Scalar Scalar;
37 typedef Matrix<Scalar,Dynamic,1> VectorType;
38 RealScalar tol = tol_error;
42 VectorType r = rhs - mat * x;
45 RealScalar r0_sqnorm = r0.squaredNorm();
46 RealScalar rhs_sqnorm = rhs.squaredNorm();
56 VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
57 VectorType y(n), z(n);
58 VectorType kt(n), ks(n);
60 VectorType s(n), t(n);
62 RealScalar tol2 = tol*tol;
63 RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
67 while ( r.squaredNorm()/rhs_sqnorm > tol2 && i<maxIters )
72 if (abs(rho) < eps2*r0_sqnorm)
77 rho = r0_sqnorm = r.squaredNorm();
81 Scalar beta = (rho/rho_old) * (alpha / w);
82 p = r + beta * (p - w * v);
86 v.noalias() = mat * y;
88 alpha = rho / r0.dot(v);
92 t.noalias() = mat * z;
94 RealScalar tmp = t.squaredNorm();
99 x += alpha * y + w * z;
103 tol_error = sqrt(r.squaredNorm()/rhs_sqnorm);
110 template<
typename _MatrixType,
111 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
116 template<
typename _MatrixType,
typename _Preconditioner>
117 struct traits<
BiCGSTAB<_MatrixType,_Preconditioner> >
119 typedef _MatrixType MatrixType;
120 typedef _Preconditioner Preconditioner;
158 template<
typename _MatrixType,
typename _Preconditioner>
159 class BiCGSTAB :
public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
161 typedef IterativeSolverBase<BiCGSTAB> Base;
162 using Base::mp_matrix;
164 using Base::m_iterations;
166 using Base::m_isInitialized;
168 typedef _MatrixType MatrixType;
169 typedef typename MatrixType::Scalar Scalar;
170 typedef typename MatrixType::Index Index;
171 typedef typename MatrixType::RealScalar RealScalar;
172 typedef _Preconditioner Preconditioner;
189 template<
typename MatrixDerived>
199 template<
typename Rhs,
typename Guess>
200 inline const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess>
203 eigen_assert(m_isInitialized &&
"BiCGSTAB is not initialized.");
204 eigen_assert(Base::rows()==b.rows()
205 &&
"BiCGSTAB::solve(): invalid number of rows of the right hand side matrix b");
206 return internal::solve_retval_with_guess
207 <
BiCGSTAB, Rhs, Guess>(*
this, b.derived(), x0);
211 template<
typename Rhs,
typename Dest>
212 void _solveWithGuess(
const Rhs& b, Dest& x)
const
215 for(
int j=0; j<b.cols(); ++j)
218 m_error = Base::m_tolerance;
220 typename Dest::ColXpr xj(x,j);
221 if(!internal::bicgstab(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
225 : m_error <= Base::m_tolerance ?
Success
227 m_isInitialized =
true;
231 template<
typename Rhs,
typename Dest>
232 void _solve(
const Rhs& b, Dest& x)
const
236 _solveWithGuess(b,x);
246 template<
typename _MatrixType,
typename _Preconditioner,
typename Rhs>
247 struct solve_retval<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
248 : solve_retval_base<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
250 typedef BiCGSTAB<_MatrixType, _Preconditioner> Dec;
251 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
253 template<typename Dest>
void evalTo(Dest& dst)
const
255 dec()._solve(rhs(),dst);
263 #endif // EIGEN_BICGSTAB_H
const internal::solve_retval_with_guess< BiCGSTAB, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: BiCGSTAB.h:201
Definition: Constants.h:378
BiCGSTAB()
Definition: BiCGSTAB.h:177
Definition: EigenBase.h:26
BiCGSTAB(const EigenBase< MatrixDerived > &A)
Definition: BiCGSTAB.h:190
Definition: Constants.h:380
Definition: Constants.h:376
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:112
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:21
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
int maxIterations() const
Definition: IterativeSolverBase.h:141