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) {
34 typedef typename Dest::RealScalar RealScalar;
35 typedef typename Dest::Scalar Scalar;
36 typedef Matrix<Scalar, Dynamic, 1> VectorType;
37 RealScalar tol = tol_error;
38 Index maxIters = iters;
41 VectorType r = rhs - mat * x;
44 RealScalar r0_norm = r0.stableNorm();
45 RealScalar r_norm = r0_norm;
46 RealScalar rhs_norm = rhs.stableNorm();
55 VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
56 VectorType y(n), z(n);
57 VectorType kt(n), ks(n);
59 VectorType s(n), t(n);
61 RealScalar eps = NumTraits<Scalar>::epsilon();
65 while (r_norm > tol && i < maxIters) {
68 if (Eigen::numext::abs(rho) / Eigen::numext::maxi(r0_norm, r_norm) < eps * Eigen::numext::mini(r0_norm, r_norm)) {
73 rho = r.squaredNorm();
74 r0_norm = r.stableNorm();
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;
85 Scalar theta = r0.dot(v);
87 RealScalar v_norm = v.stableNorm();
88 if (Eigen::numext::abs(theta) / Eigen::numext::maxi(r0_norm, v_norm) < eps * Eigen::numext::mini(r0_norm, v_norm)) {
91 r0_norm = r0.stableNorm();
95 if (restarts++ == 0) i = 0;
101 z = precond.solve(s);
102 t.noalias() = mat * z;
104 RealScalar tmp = t.squaredNorm();
105 if (tmp > RealScalar(0)) {
110 x += alpha * y + w * z;
112 r_norm = r.stableNorm();
116 tol_error = r_norm / rhs_norm;
123template <
typename MatrixType_,
typename Preconditioner_ = DiagonalPreconditioner<
typename MatrixType_::Scalar> >
128template <
typename MatrixType_,
typename Preconditioner_>
129struct traits<BiCGSTAB<MatrixType_, Preconditioner_> > {
130 typedef MatrixType_ MatrixType;
131 typedef Preconditioner_ Preconditioner;
167template <
typename MatrixType_,
typename Preconditioner_>
172 using Base::m_isInitialized;
173 using Base::m_iterations;
177 typedef MatrixType_ MatrixType;
178 typedef typename MatrixType::Scalar Scalar;
179 typedef typename MatrixType::RealScalar RealScalar;
180 typedef Preconditioner_ Preconditioner;
196 template <
typename MatrixDerived>
202 template <
typename Rhs,
typename Dest>
203 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
205 m_error = Base::m_tolerance;
207 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:168
BiCGSTAB(const EigenBase< MatrixDerived > &A)
Definition BiCGSTAB.h:197
BiCGSTAB()
Definition BiCGSTAB.h:184
IterativeSolverBase()
Definition IterativeSolverBase.h:142
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 B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Definition EigenBase.h:33