// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2011-2014 Gael Guennebaud // Copyright (C) 2012 Désiré Nuentsa-Wakam // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_BICGSTAB_H #define EIGEN_BICGSTAB_H // IWYU pragma: private #include "./InternalHeaderCheck.h" namespace Eigen { namespace internal { /** \internal Low-level bi conjugate gradient stabilized algorithm * \param mat The matrix A * \param rhs The right hand side vector b * \param x On input and initial solution, on output the computed solution. * \param precond A preconditioner being able to efficiently solve for an * approximation of Ax=b (regardless of b) * \param iters On input the max number of iteration, on output the number of performed iterations. * \param tol_error On input the tolerance error, on output an estimation of the relative error. * \return false in the case of numerical issue, for example a break down of BiCGSTAB. */ template bool bicgstab(const MatrixType& mat, const Rhs& rhs, Dest& x, const Preconditioner& precond, Index& iters, typename Dest::RealScalar& tol_error) { typedef typename Dest::RealScalar RealScalar; typedef typename Dest::Scalar Scalar; typedef Matrix VectorType; RealScalar tol = tol_error; Index maxIters = iters; Index n = mat.cols(); VectorType r = rhs - mat * x; VectorType r0 = r; RealScalar r0_norm = r0.stableNorm(); RealScalar r_norm = r0_norm; RealScalar rhs_norm = rhs.stableNorm(); if (rhs_norm == 0) { x.setZero(); return true; } Scalar rho(1); Scalar alpha(0); Scalar w(1); VectorType v = VectorType::Zero(n), p = VectorType::Zero(n); VectorType y(n), z(n); VectorType kt(n), ks(n); VectorType s(n), t(n); RealScalar eps = NumTraits::epsilon(); Index i = 0; Index restarts = 0; while (r_norm > tol && i < maxIters) { Scalar rho_old = rho; rho = r0.dot(r); if (Eigen::numext::abs(rho) / Eigen::numext::maxi(r0_norm, r_norm) < eps * Eigen::numext::mini(r0_norm, r_norm)) { // The new residual vector became too orthogonal to the arbitrarily chosen direction r0 // Let's restart with a new r0: r = rhs - mat * x; r0 = r; rho = r.squaredNorm(); r0_norm = r.stableNorm(); alpha = Scalar(0); w = Scalar(1); if (restarts++ == 0) i = 0; } Scalar beta = (rho / rho_old) * (alpha / w); p = r + beta * (p - w * v); y = precond.solve(p); v.noalias() = mat * y; Scalar theta = r0.dot(v); // For small angles ∠(r0, v) < eps, random restart. RealScalar v_norm = v.stableNorm(); if (Eigen::numext::abs(theta) / Eigen::numext::maxi(r0_norm, v_norm) < eps * Eigen::numext::mini(r0_norm, v_norm)) { r = rhs - mat * x; r0.setRandom(); r0_norm = r0.stableNorm(); rho = Scalar(1); alpha = Scalar(0); w = Scalar(1); if (restarts++ == 0) i = 0; continue; } alpha = rho / theta; s = r - alpha * v; z = precond.solve(s); t.noalias() = mat * z; RealScalar tmp = t.squaredNorm(); if (tmp > RealScalar(0)) { w = t.dot(s) / tmp; } else { w = Scalar(0); } x += alpha * y + w * z; r = s - w * t; r_norm = r.stableNorm(); ++i; } tol_error = r_norm / rhs_norm; iters = i; return true; } } // namespace internal template > class BiCGSTAB; namespace internal { template struct traits > { typedef MatrixType_ MatrixType; typedef Preconditioner_ Preconditioner; }; } // namespace internal /** \ingroup IterativeLinearSolvers_Module * \brief A bi conjugate gradient stabilized solver for sparse square problems * * This class allows to solve for A.x = b sparse linear problems using a bi conjugate gradient * stabilized algorithm. The vectors x and b can be either dense or sparse. * * \tparam MatrixType_ the type of the sparse matrix A, can be a dense or a sparse matrix. * \tparam Preconditioner_ the type of the preconditioner. Default is DiagonalPreconditioner * * \implsparsesolverconcept * * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations * and NumTraits::epsilon() for the tolerance. * * The tolerance corresponds to the relative residual error: |Ax-b|/|b| * * \b Performance: when using sparse matrices, best performance is achied for a row-major sparse matrix format. * Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled. * See \ref TopicMultiThreading for details. * * This class can be used as the direct solver classes. Here is a typical usage example: * \include BiCGSTAB_simple.cpp * * By default the iterations start with x=0 as an initial guess of the solution. * One can control the start using the solveWithGuess() method. * * BiCGSTAB can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink. * * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner */ template class BiCGSTAB : public IterativeSolverBase > { typedef IterativeSolverBase Base; using Base::m_error; using Base::m_info; using Base::m_isInitialized; using Base::m_iterations; using Base::matrix; public: typedef MatrixType_ MatrixType; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; typedef Preconditioner_ Preconditioner; public: /** Default constructor. */ BiCGSTAB() : Base() {} /** Initialize the solver with matrix \a A for further \c Ax=b solving. * * This constructor is a shortcut for the default constructor followed * by a call to compute(). * * \warning this class stores a reference to the matrix A as well as some * precomputed values that depend on it. Therefore, if \a A is changed * this class becomes invalid. Call compute() to update it with the new * matrix A, or modify a copy of A. */ template explicit BiCGSTAB(const EigenBase& A) : Base(A.derived()) {} ~BiCGSTAB() {} /** \internal */ template void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const { m_iterations = Base::maxIterations(); m_error = Base::m_tolerance; bool ret = internal::bicgstab(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error); m_info = (!ret) ? NumericalIssue : m_error <= Base::m_tolerance ? Success : NoConvergence; } protected: }; } // end namespace Eigen #endif // EIGEN_BICGSTAB_H