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;
    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;
    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,
   117 template< 
typename _MatrixType, 
typename _Preconditioner>
   157 template< 
typename _MatrixType, 
typename _Preconditioner>
   163   using Base::m_iterations;
   165   using Base::m_isInitialized;
   168   typedef typename MatrixType::Scalar 
Scalar;
   187   template<
typename MatrixDerived>
   193   template<
typename Rhs,
typename Dest>
   197     for(
Index j=0; j<b.cols(); ++j)
   199       m_iterations = Base::maxIterations();
   200       m_error = Base::m_tolerance;
   203       if(!
internal::bicgstab(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error))
   207            : m_error <= Base::m_tolerance ? 
Success   209     m_isInitialized = 
true;
   213   using Base::_solve_impl;
   214   template<
typename Rhs,
typename Dest>
   217     x.resize(this->rows(),b.cols());
   219     _solve_with_guess_impl(b,x);
   228 #endif // EIGEN_BICGSTAB_H A preconditioner based on the digonal entries. 
bool bicgstab(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
void _solve_with_guess_impl(const Rhs &b, Dest &x) const
_Preconditioner Preconditioner
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1, !IsRowMajor > ColXpr
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
void _solve_impl(const MatrixBase< Rhs > &b, Dest &x) const
double rho(Eigen::Vector2d xy)
BiCGSTAB(const EigenBase< MatrixDerived > &A)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API. 
_Preconditioner Preconditioner
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index cols() const
MatrixType::RealScalar RealScalar
A bi conjugate gradient stabilized solver for sparse square problems. 
MatrixType::Scalar Scalar
IterativeSolverBase< BiCGSTAB > Base
Base class for linear iterative solvers. 
Base class for all dense matrices, vectors, and expressions.