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, Index& 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;
39 Index maxIters = iters;
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*rhs_sqnorm;
63 RealScalar eps2 = NumTraits<Scalar>::epsilon()*NumTraits<Scalar>::epsilon();
67 while ( r.squaredNorm() > tol2 && i<maxIters )
72 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;
153 template<
typename _MatrixType,
typename _Preconditioner>
154 class BiCGSTAB :
public IterativeSolverBase<BiCGSTAB<_MatrixType,_Preconditioner> >
156 typedef IterativeSolverBase<BiCGSTAB> Base;
157 using Base::mp_matrix;
159 using Base::m_iterations;
161 using Base::m_isInitialized;
163 typedef _MatrixType MatrixType;
164 typedef typename MatrixType::Scalar Scalar;
165 typedef typename MatrixType::RealScalar RealScalar;
166 typedef _Preconditioner Preconditioner;
183 explicit BiCGSTAB(
const MatrixType& A) : Base(A) {}
188 template<
typename Rhs,
typename Dest>
189 void _solve_with_guess_impl(
const Rhs& b, Dest& x)
const
192 for(Index j=0; j<b.cols(); ++j)
195 m_error = Base::m_tolerance;
197 typename Dest::ColXpr xj(x,j);
198 if(!internal::bicgstab(mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
202 : m_error <= Base::m_tolerance ?
Success
204 m_isInitialized =
true;
208 using Base::_solve_impl;
209 template<
typename Rhs,
typename Dest>
210 void _solve_impl(
const MatrixBase<Rhs>& b, Dest& x)
const
212 x.resize(this->rows(),b.cols());
214 _solve_with_guess_impl(b,x);
223 #endif // EIGEN_BICGSTAB_H
BiCGSTAB(const MatrixType &A)
Definition: BiCGSTAB.h:183
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