11#ifndef EIGEN_BICGSTAB_H
12#define EIGEN_BICGSTAB_H
28template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
29bool 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);
111template<
typename _MatrixType,
112 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
117template<
typename _MatrixType,
typename _Preconditioner>
118struct traits<BiCGSTAB<_MatrixType,_Preconditioner> >
120 typedef _MatrixType MatrixType;
121 typedef _Preconditioner Preconditioner;
157template<
typename _MatrixType,
typename _Preconditioner>
163 using Base::m_iterations;
165 using Base::m_isInitialized;
167 typedef _MatrixType MatrixType;
168 typedef typename MatrixType::Scalar Scalar;
169 typedef typename MatrixType::RealScalar RealScalar;
170 typedef _Preconditioner Preconditioner;
187 template<
typename MatrixDerived>
193 template<
typename Rhs,
typename Dest>
194 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
197 m_error = Base::m_tolerance;
199 bool ret = internal::bicgstab(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
202 : m_error <= Base::m_tolerance ?
Success
A bi conjugate gradient stabilized solver for sparse square problems.
Definition: BiCGSTAB.h:159
BiCGSTAB(const EigenBase< MatrixDerived > &A)
Definition: BiCGSTAB.h:188
BiCGSTAB()
Definition: BiCGSTAB.h:175
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:144
Index maxIterations() const
Definition: IterativeSolverBase.h:281
@ NumericalIssue
Definition: Constants.h:444
@ Success
Definition: Constants.h:442
@ NoConvergence
Definition: Constants.h:446
Namespace containing all symbols from the Eigen library.
Definition: Core:141
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Definition: EigenBase.h:30