11#ifndef EIGEN_BICGSTAB_H
12#define EIGEN_BICGSTAB_H
15#include "./InternalHeaderCheck.h"
31template <
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
32bool bicgstab(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond, Index& iters,
33 typename Dest::RealScalar& tol_error) {
36 typedef typename Dest::RealScalar RealScalar;
37 typedef typename Dest::Scalar Scalar;
38 typedef Matrix<Scalar, Dynamic, 1> VectorType;
39 RealScalar tol = tol_error;
40 Index maxIters = iters;
43 VectorType r = rhs - mat * x;
46 RealScalar r0_sqnorm = r0.squaredNorm();
47 RealScalar rhs_sqnorm = rhs.squaredNorm();
48 if (rhs_sqnorm == 0) {
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) {
71 if (
abs(rho) < eps2 * r0_sqnorm) {
76 rho = r0_sqnorm = r.squaredNorm();
77 if (restarts++ == 0) i = 0;
79 Scalar beta = (rho / rho_old) * (alpha / w);
80 p = r + beta * (p - w * v);
84 v.noalias() = mat * y;
86 alpha = rho / r0.dot(v);
90 t.noalias() = mat * z;
92 RealScalar tmp = t.squaredNorm();
93 if (tmp > RealScalar(0))
97 x += alpha * y + w * z;
101 tol_error =
sqrt(r.squaredNorm() / rhs_sqnorm);
108template <
typename MatrixType_,
typename Preconditioner_ = DiagonalPreconditioner<
typename MatrixType_::Scalar> >
113template <
typename MatrixType_,
typename Preconditioner_>
114struct traits<BiCGSTAB<MatrixType_, Preconditioner_> > {
115 typedef MatrixType_ MatrixType;
116 typedef Preconditioner_ Preconditioner;
152template <
typename MatrixType_,
typename Preconditioner_>
157 using Base::m_isInitialized;
158 using Base::m_iterations;
162 typedef MatrixType_ MatrixType;
163 typedef typename MatrixType::Scalar Scalar;
164 typedef typename MatrixType::RealScalar RealScalar;
165 typedef Preconditioner_ Preconditioner;
181 template <
typename MatrixDerived>
187 template <
typename Rhs,
typename Dest>
188 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
190 m_error = Base::m_tolerance;
192 bool ret = internal::bicgstab(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
A bi conjugate gradient stabilized solver for sparse square problems.
Definition BiCGSTAB.h:153
BiCGSTAB(const EigenBase< MatrixDerived > &A)
Definition BiCGSTAB.h:182
BiCGSTAB()
Definition BiCGSTAB.h:169
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:124
Index maxIterations() const
Definition IterativeSolverBase.h:251
@ NumericalIssue
Definition Constants.h:442
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33