bicgstabsolver.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_BICG_STAB_SOLVER_HH
28#define EWOMS_BICG_STAB_SOLVER_HH
29
32#include "linearsolverreport.hh"
33
36
37#include <opm/common/Exceptions.hpp>
38
39#include <memory>
40
41namespace Opm {
42namespace Linear {
52template <class LinearOperator, class Vector, class Preconditioner>
54{
56 using Scalar = typename LinearOperator::field_type;
57
58public:
59 BiCGStabSolver(Preconditioner& preconditioner,
60 ConvergenceCriterion& convergenceCriterion,
61 Dune::ScalarProduct<Vector>& scalarProduct)
62 : preconditioner_(preconditioner)
63 , convergenceCriterion_(convergenceCriterion)
64 , scalarProduct_(scalarProduct)
65 {
66 A_ = nullptr;
67 b_ = nullptr;
68
69 maxIterations_ = 1000;
70 }
71
76 void setMaxIterations(unsigned value)
77 { maxIterations_ = value; }
78
83 unsigned maxIterations() const
84 { return maxIterations_; }
85
96 void setVerbosity(unsigned value)
97 { verbosity_ = value; }
98
102 unsigned verbosity() const
103 { return verbosity_; }
104
108 void setLinearOperator(const LinearOperator* A)
109 { A_ = A; }
110
114 void setRhs(const Vector* b)
115 { b_ = b; }
116
120 bool apply(Vector& x)
121 {
122 // epsilon used for detecting breakdowns
123 const Scalar breakdownEps = std::numeric_limits<Scalar>::min() * Scalar(1e10);
124
125 // start the stop watch for the solution proceedure, but make sure that it is
126 // turned off regardless of how we leave the stadium. (i.e., that the timer gets
127 // stopped in case exceptions are thrown as well as if the method returns
128 // regularly.)
129 report_.reset();
130 TimerGuard reportTimerGuard(report_.timer());
131 report_.timer().start();
132
133 // preconditioned stabilized biconjugate gradient method
134 //
135 // See https://en.wikipedia.org/wiki/Biconjugate_gradient_stabilized_method,
136 // (article date: December 19, 2016)
137
138 // set the initial solution to the zero vector
139 x = 0.0;
140
141 // prepare the preconditioner. to allow some optimizations, we assume that the
142 // preconditioner does not change the initial solution x if the initial solution
143 // is a zero vector.
144 Vector r = *b_;
145 preconditioner_.pre(x, r);
146
147#ifndef NDEBUG
148 // ensure that the preconditioner does not change the initial solution. since
149 // this is a debugging check, we don't care if it does not work properly in
150 // parallel. (because this goes wrong, it should be considered to be a bug in the
151 // code anyway.)
152 for (unsigned i = 0; i < x.size(); ++i) {
153 const auto& u = x[i];
154 if (u*u != 0.0)
155 throw std::logic_error("The preconditioner is assumed not to modify the initial solution!");
156 }
157#endif // NDEBUG
158
159 convergenceCriterion_.setInitial(x, r);
160 if (convergenceCriterion_.converged()) {
161 report_.setConverged(true);
162 return report_.converged();
163 }
164
165 if (verbosity_ > 0) {
166 std::cout << "-------- BiCGStabSolver --------" << std::endl;
167 convergenceCriterion_.printInitial();
168 }
169
170 // r0 = b - Ax (i.e., r -= A*x_0 = b, because x_0 == 0)
171 //A_->applyscaleadd(/*alpha=*/-1.0, x, r);
172
173 // r0hat = r0
174 const Vector& r0hat = *b_;
175
176 // rho0 = alpha = omega0 = 1
177 Scalar rho = 1.0;
178 Scalar alpha = 1.0;
179 Scalar omega = 1.0;
180
181 // v_0 = p_0 = 0;
182 Vector v(r);
183 v = 0.0;
184 Vector p(v);
185
186 // create all the temporary vectors which we need. Be aware that some of them
187 // actually point to the same object because they are not needed at the same time!
188 Vector y(x);
189 Vector& h(x);
190 Vector& s(r);
191 Vector z(x);
192 Vector& t(y);
193 unsigned n = x.size();
194
195 for (; report_.iterations() < maxIterations_; report_.increment()) {
196 // rho_i = (r0hat,r_(i-1))
197 Scalar rho_i = scalarProduct_.dot(r0hat, r);
198
199 // beta = (rho_i/rho_(i-1))*(alpha/omega_(i-1))
200 if (std::abs(rho) <= breakdownEps || std::abs(omega) <= breakdownEps)
201 throw NumericalProblem("Breakdown of the BiCGStab solver (division by zero)");
202 Scalar beta = (rho_i/rho)*(alpha/omega);
203
204 // make rho correspond to the current iteration (i.e., forget rho_(i-1))
205 rho = rho_i;
206
207 // this loop conflates the following operations:
208 //
209 // p_i = r_(i-1) + beta*(p_(i-1) - omega_(i-1)*v_(i-1))
210 // y = p
211 for (unsigned i = 0; i < n; ++i) {
212 // p_i = r_(i-1) + beta*(p_(i-1) - omega_(i-1)*v_(i-1))
213 auto tmp = v[i];
214 tmp *= omega;
215 tmp -= p[i];
216 tmp *= -beta;
217 p[i] = r[i];
218 p[i] += tmp;
219
220 // y = p; not required because the precontioner overwrites y anyway...
221 // y[i] = p[i];
222 }
223
224 // y = K^-1 * p_i
225 preconditioner_.apply(y, p);
226
227 // v_i = A*y
228 A_->apply(y, v);
229
230 // alpha = rho_i/(r0hat,v_i)
231 Scalar denom = scalarProduct_.dot(r0hat, v);
232 if (std::abs(denom) <= breakdownEps)
233 throw NumericalProblem("Breakdown of the BiCGStab solver (division by zero)");
234 alpha = rho_i/denom;
235 if (std::abs(alpha) <= breakdownEps)
236 throw NumericalProblem("Breakdown of the BiCGStab solver (stagnation detected)");
237
238 // h = x_(i-1) + alpha*y
239 // s = r_(i-1) - alpha*v_i
240 for (unsigned i = 0; i < n; ++i) {
241 auto tmp = y[i];
242 tmp *= alpha;
243 tmp += x[i];
244 h[i] = tmp;
245
246 //s[i] = r[i]; // not necessary because r and s are the same object
247 tmp = v[i];
248 tmp *= alpha;
249 s[i] -= tmp;
250 }
251
252 // do convergence check and print terminal output
253 convergenceCriterion_.update(/*curSol=*/h, /*delta=*/y, s);
254 if (convergenceCriterion_.converged()) {
255 if (verbosity_ > 0) {
256 convergenceCriterion_.print(report_.iterations() + 0.5);
257 std::cout << "-------- /BiCGStabSolver --------" << std::endl;
258 }
259
260 // x = h; // not necessary because x and h are the same object
261 preconditioner_.post(x);
262 report_.setConverged(true);
263 return report_.converged();
264 }
265 else if (convergenceCriterion_.failed()) {
266 if (verbosity_ > 0) {
267 convergenceCriterion_.print(report_.iterations() + 0.5);
268 std::cout << "-------- /BiCGStabSolver --------" << std::endl;
269 }
270
271 report_.setConverged(false);
272 return report_.converged();
273 }
274
275 if (verbosity_ > 1)
276 convergenceCriterion_.print(report_.iterations() + 0.5);
277
278 // z = K^-1*s
279 z = s;
280 preconditioner_.apply(z, s);
281
282 // t = Az
283 t = z;
284 A_->apply(z, t);
285
286 // omega_i = (t*s)/(t*t)
287 denom = scalarProduct_.dot(t, t);
288 if (std::abs(denom) <= breakdownEps)
289 throw NumericalProblem("Breakdown of the BiCGStab solver (division by zero)");
290 omega = scalarProduct_.dot(t, s)/denom;
291 if (std::abs(omega) <= breakdownEps)
292 throw NumericalProblem("Breakdown of the BiCGStab solver (stagnation detected)");
293
294 // x_i = h + omega_i*z
295 // x = h; // not necessary because x and h are the same object
296 x.axpy(/*a=*/omega, /*y=*/z);
297
298 // do convergence check and print terminal output
299 convergenceCriterion_.update(/*curSol=*/x, /*delta=*/z, r);
300 if (convergenceCriterion_.converged()) {
301 if (verbosity_ > 0) {
302 convergenceCriterion_.print(1.0 + report_.iterations());
303 std::cout << "-------- /BiCGStabSolver --------" << std::endl;
304 }
305
306 preconditioner_.post(x);
307 report_.setConverged(true);
308 return report_.converged();
309 }
310 else if (convergenceCriterion_.failed()) {
311 if (verbosity_ > 0) {
312 convergenceCriterion_.print(1.0 + report_.iterations());
313 std::cout << "-------- /BiCGStabSolver --------" << std::endl;
314 }
315
316 report_.setConverged(false);
317 return report_.converged();
318 }
319
320 if (verbosity_ > 1)
321 convergenceCriterion_.print(1.0 + report_.iterations());
322
323 // r_i = s - omega*t
324 // r = s; // not necessary because r and s are the same object
325 r.axpy(/*a=*/-omega, /*y=*/t);
326 }
327
328 report_.setConverged(false);
329 return report_.converged();
330 }
331
333 {
334 convergenceCriterion_ = &crit;
335 }
336
337 const SolverReport& report() const
338 { return report_; }
339
340private:
341 const LinearOperator* A_;
342 const Vector* b_;
343
344 Preconditioner& preconditioner_;
345 ConvergenceCriterion& convergenceCriterion_;
346 Dune::ScalarProduct<Vector>& scalarProduct_;
347 SolverReport report_;
348
349 unsigned maxIterations_;
350 unsigned verbosity_;
351};
352
353} // namespace Linear
354} // namespace Opm
355
356#endif
Implements a preconditioned stabilized BiCG linear solver.
Definition: bicgstabsolver.hh:54
void setConvergenceCriterion(ConvergenceCriterion &crit)
Definition: bicgstabsolver.hh:332
void setVerbosity(unsigned value)
Set the verbosity level of the linear solver.
Definition: bicgstabsolver.hh:96
void setLinearOperator(const LinearOperator *A)
Set the matrix "A" of the linear system.
Definition: bicgstabsolver.hh:108
bool apply(Vector &x)
Run the stabilized BiCG solver and store the result into the "x" vector.
Definition: bicgstabsolver.hh:120
unsigned verbosity() const
Return the verbosity level of the linear solver.
Definition: bicgstabsolver.hh:102
void setMaxIterations(unsigned value)
Set the maximum number of iterations before we give up without achieving convergence.
Definition: bicgstabsolver.hh:76
unsigned maxIterations() const
Return the maximum number of iterations before we give up without achieving convergence.
Definition: bicgstabsolver.hh:83
const SolverReport & report() const
Definition: bicgstabsolver.hh:337
void setRhs(const Vector *b)
Set the right hand side "b" of the linear system.
Definition: bicgstabsolver.hh:114
BiCGStabSolver(Preconditioner &preconditioner, ConvergenceCriterion &convergenceCriterion, Dune::ScalarProduct< Vector > &scalarProduct)
Definition: bicgstabsolver.hh:59
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:56
virtual void setInitial(const Vector &curSol, const Vector &curResid)=0
Set the initial solution of the linear system of equations.
virtual void update(const Vector &curSol, const Vector &changeIndicator, const Vector &curResid)=0
Update the internal members of the convergence criterion with the current solution.
virtual bool failed() const
Returns true if the convergence criterion cannot be met anymore because the solver has broken down.
Definition: convergencecriterion.hh:115
virtual bool converged() const =0
Returns true if and only if the convergence criterion is met.
virtual void printInitial(std::ostream &=std::cout) const
Prints the initial information about the convergence behaviour.
Definition: convergencecriterion.hh:135
virtual void print(Scalar, std::ostream &=std::cout) const
Prints the information about the convergence behaviour for the current iteration.
Definition: convergencecriterion.hh:146
Collects summary information about the execution of the linear solver.
Definition: linearsolverreport.hh:42
bool converged() const
Definition: linearsolverreport.hh:69
unsigned iterations() const
Definition: linearsolverreport.hh:60
void reset()
Definition: linearsolverreport.hh:47
void setConverged(bool value)
Definition: linearsolverreport.hh:72
const Opm::Timer & timer() const
Definition: linearsolverreport.hh:54
void increment()
Definition: linearsolverreport.hh:63
A simple class which makes sure that a timer gets stopped if an exception is thrown.
Definition: timerguard.hh:41
void start()
Start counting the time resources used by the simulation.
Definition: timer.hh:62
Define some base class for the convergence criteria of the linear solvers of DUNE-ISTL.
Definition: blackoilboundaryratevector.hh:37