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;
43 VectorType r = rhs - mat * x;
46 RealScalar r0_sqnorm = r0.squaredNorm();
47 RealScalar rhs_sqnorm = rhs.squaredNorm();
57 VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
58 VectorType y(n), z(n);
59 VectorType kt(n), ks(n);
61 VectorType s(n), t(n);
63 RealScalar tol2 = tol*tol;
64 RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
68 while ( r.squaredNorm()/rhs_sqnorm > tol2 && i<maxIters )
73 if (abs(rho) < eps2*r0_sqnorm)
78 rho = r0_sqnorm = r.squaredNorm();
82 Scalar beta = (rho/rho_old) * (alpha / w);
83 p = r + beta * (p - w * v);
87 v.noalias() = mat * y;
89 alpha = rho / r0.dot(v);
93 t.noalias() = mat * z;
95 RealScalar tmp = t.squaredNorm();
100 x += alpha * y + w * z;
104 tol_error = sqrt(r.squaredNorm()/rhs_sqnorm);
111 template<
typename _MatrixType,
112 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
117 template<
typename _MatrixType,
typename _Preconditioner>
118 struct traits<
BiCGSTAB<_MatrixType,_Preconditioner> >
120 typedef _MatrixType MatrixType;
121 typedef _Preconditioner Preconditioner;
172 template<
typename _MatrixType,
typename _Preconditioner>
173 class BiCGSTAB :
public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
175 typedef IterativeSolverBase<BiCGSTAB> Base;
176 using Base::mp_matrix;
178 using Base::m_iterations;
180 using Base::m_isInitialized;
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;
212 template<
typename Rhs,
typename Guess>
213 inline const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess>
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);
224 template<
typename Rhs,
typename Dest>
225 void _solveWithGuess(
const Rhs& b, Dest& x)
const
228 for(
int j=0; j<b.cols(); ++j)
231 m_error = Base::m_tolerance;
233 typename Dest::ColXpr xj(x,j);
234 if(!internal::bicgstab(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
238 : m_error <= Base::m_tolerance ?
Success
240 m_isInitialized =
true;
244 template<
typename Rhs,
typename Dest>
245 void _solve(
const Rhs& b, Dest& x)
const
249 _solveWithGuess(b,x);
259 template<
typename _MatrixType,
typename _Preconditioner,
typename Rhs>
260 struct solve_retval<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
261 : solve_retval_base<BiCGSTAB<_MatrixType, _Preconditioner>, Rhs>
263 typedef BiCGSTAB<_MatrixType, _Preconditioner> Dec;
264 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
266 template<typename Dest>
void evalTo(Dest& dst)
const
268 dec()._solve(rhs(),dst);
276 #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:214
Definition: Constants.h:378
BiCGSTAB()
Definition: BiCGSTAB.h:191
BiCGSTAB(const MatrixType &A)
Definition: BiCGSTAB.h:203
Definition: Constants.h:380
Definition: Constants.h:376
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:113
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:136