Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
BiCGSTAB.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_BICGSTAB_H
12#define EIGEN_BICGSTAB_H
13
14// IWYU pragma: private
15#include "./InternalHeaderCheck.h"
16
17namespace Eigen {
18
19namespace internal {
20
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;
39
40 Index n = mat.cols();
41 VectorType r = rhs - mat * x;
42 VectorType r0 = r;
43
44 RealScalar r0_norm = r0.stableNorm();
45 RealScalar r_norm = r0_norm;
46 RealScalar rhs_norm = rhs.stableNorm();
47 if (rhs_norm == 0) {
48 x.setZero();
49 return true;
50 }
51 Scalar rho(1);
52 Scalar alpha(0);
53 Scalar w(1);
54
55 VectorType v = VectorType::Zero(n), p = VectorType::Zero(n);
56 VectorType y(n), z(n);
57 VectorType kt(n), ks(n);
58
59 VectorType s(n), t(n);
60
61 RealScalar eps = NumTraits<Scalar>::epsilon();
62 Index i = 0;
63 Index restarts = 0;
64
65 while (r_norm > tol && i < maxIters) {
66 Scalar rho_old = rho;
67 rho = r0.dot(r);
68 if (Eigen::numext::abs(rho) / Eigen::numext::maxi(r0_norm, r_norm) < eps * Eigen::numext::mini(r0_norm, r_norm)) {
69 // The new residual vector became too orthogonal to the arbitrarily chosen direction r0
70 // Let's restart with a new r0:
71 r = rhs - mat * x;
72 r0 = r;
73 rho = r.squaredNorm();
74 r0_norm = r.stableNorm();
75 alpha = Scalar(0);
76 w = Scalar(1);
77 if (restarts++ == 0) i = 0;
78 }
79 Scalar beta = (rho / rho_old) * (alpha / w);
80 p = r + beta * (p - w * v);
81
82 y = precond.solve(p);
83
84 v.noalias() = mat * y;
85 Scalar theta = r0.dot(v);
86 // For small angles ∠(r0, v) < eps, random restart.
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)) {
89 r = rhs - mat * x;
90 r0.setRandom();
91 r0_norm = r0.stableNorm();
92 rho = Scalar(1);
93 alpha = Scalar(0);
94 w = Scalar(1);
95 if (restarts++ == 0) i = 0;
96 continue;
97 }
98 alpha = rho / theta;
99 s = r - alpha * v;
100
101 z = precond.solve(s);
102 t.noalias() = mat * z;
103
104 RealScalar tmp = t.squaredNorm();
105 if (tmp > RealScalar(0)) {
106 w = t.dot(s) / tmp;
107 } else {
108 w = Scalar(0);
109 }
110 x += alpha * y + w * z;
111 r = s - w * t;
112 r_norm = r.stableNorm();
113 ++i;
114 }
115
116 tol_error = r_norm / rhs_norm;
117 iters = i;
118 return true;
119}
120
121} // namespace internal
122
123template <typename MatrixType_, typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
124class BiCGSTAB;
125
126namespace internal {
127
128template <typename MatrixType_, typename Preconditioner_>
129struct traits<BiCGSTAB<MatrixType_, Preconditioner_> > {
130 typedef MatrixType_ MatrixType;
131 typedef Preconditioner_ Preconditioner;
132};
133
134} // namespace internal
135
167template <typename MatrixType_, typename Preconditioner_>
168class BiCGSTAB : public IterativeSolverBase<BiCGSTAB<MatrixType_, Preconditioner_> > {
170 using Base::m_error;
171 using Base::m_info;
172 using Base::m_isInitialized;
173 using Base::m_iterations;
174 using Base::matrix;
175
176 public:
177 typedef MatrixType_ MatrixType;
178 typedef typename MatrixType::Scalar Scalar;
179 typedef typename MatrixType::RealScalar RealScalar;
180 typedef Preconditioner_ Preconditioner;
181
182 public:
184 BiCGSTAB() : Base() {}
185
196 template <typename MatrixDerived>
197 explicit BiCGSTAB(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
198
199 ~BiCGSTAB() {}
200
202 template <typename Rhs, typename Dest>
203 void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const {
204 m_iterations = Base::maxIterations();
205 m_error = Base::m_tolerance;
206
207 bool ret = internal::bicgstab(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
208
209 m_info = (!ret) ? NumericalIssue : m_error <= Base::m_tolerance ? Success : NoConvergence;
210 }
211
212 protected:
213};
214
215} // end namespace Eigen
216
217#endif // EIGEN_BICGSTAB_H
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
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