Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
solver_bicgstab.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_solver_bicgstab_h
17#define dealii_solver_bicgstab_h
18
19
20#include <deal.II/base/config.h>
21
25
26#include <deal.II/lac/solver.h>
28
29#include <cmath>
30
32
35
36namespace internal
37{
43 {
44 protected:
48 double alpha;
52 double beta;
56 double omega;
60 double rho;
64 double rhobar;
65
69 unsigned int step;
70
74 double res;
75
81 };
82} // namespace internal
83
122template <typename VectorType = Vector<double>>
123class SolverBicgstab : public SolverBase<VectorType>,
125{
126public:
138 {
147 const bool exact_residual = true,
148 const double breakdown =
152 {}
160 double breakdown;
161 };
162
168 const AdditionalData & data = AdditionalData());
169
175 const AdditionalData &data = AdditionalData());
176
180 virtual ~SolverBicgstab() override = default;
181
185 template <typename MatrixType, typename PreconditionerType>
186 void
187 solve(const MatrixType & A,
188 VectorType & x,
189 const VectorType & b,
190 const PreconditionerType &preconditioner);
191
192protected:
196 VectorType *Vx;
197
202
207
212
217
222
227
232
236 const VectorType *Vb;
237
241 template <typename MatrixType>
242 double
243 criterion(const MatrixType &A, const VectorType &x, const VectorType &b);
244
250 virtual void
251 print_vectors(const unsigned int step,
252 const VectorType & x,
253 const VectorType & r,
254 const VectorType & d) const;
255
260
261private:
267 {
270 unsigned int last_step;
272
275 const unsigned int last_step,
276 const double last_residual);
277 };
278
283 template <typename MatrixType, typename PreconditionerType>
285 iterate(const MatrixType &A, const PreconditionerType &preconditioner);
286};
287
289/*-------------------------Inline functions -------------------------------*/
290
291#ifndef DOXYGEN
292
293
294template <typename VectorType>
296 const bool breakdown,
297 const SolverControl::State state,
298 const unsigned int last_step,
299 const double last_residual)
300 : breakdown(breakdown)
301 , state(state)
302 , last_step(last_step)
303 , last_residual(last_residual)
304{}
305
306
307
308template <typename VectorType>
311 const AdditionalData & data)
312 : SolverBase<VectorType>(cn, mem)
313 , Vx(nullptr)
314 , Vb(nullptr)
315 , additional_data(data)
316{}
317
318
319
320template <typename VectorType>
322 const AdditionalData &data)
323 : SolverBase<VectorType>(cn)
324 , Vx(nullptr)
325 , Vb(nullptr)
326 , additional_data(data)
327{}
328
329
330
331template <typename VectorType>
332template <typename MatrixType>
333double
335 const VectorType &x,
336 const VectorType &b)
337{
338 A.vmult(*Vt, x);
339 Vt->add(-1., b);
340 res = Vt->l2_norm();
341
342 return res;
343}
344
345
346
347template <typename VectorType>
348void
350 const VectorType &,
351 const VectorType &,
352 const VectorType &) const
353{}
354
355
356
357template <typename VectorType>
358template <typename MatrixType, typename PreconditionerType>
360SolverBicgstab<VectorType>::iterate(const MatrixType & A,
361 const PreconditionerType &preconditioner)
362{
363 A.vmult(*Vr, *Vx);
364 Vr->sadd(-1., 1., *Vb);
365 res = Vr->l2_norm();
366
367 SolverControl::State state = this->iteration_status(step, res, *Vx);
368 if (state == SolverControl::State::success)
369 return IterationResult(false, state, step, res);
370
371 alpha = omega = rho = 1.;
372
373 VectorType &r = *Vr;
374 VectorType &rbar = *Vrbar;
375 VectorType &p = *Vp;
376 VectorType &y = *Vy;
377 VectorType &z = *Vz;
378 VectorType &t = *Vt;
379 VectorType &v = *Vv;
380
381 rbar = r;
382 bool startup = true;
383
384 do
385 {
386 ++step;
387
388 rhobar = r * rbar;
389 if (std::fabs(rhobar) < additional_data.breakdown)
390 {
391 return IterationResult(true, state, step, res);
392 }
393 beta = rhobar * alpha / (rho * omega);
394 rho = rhobar;
395 if (startup == true)
396 {
397 p = r;
398 startup = false;
399 }
400 else
401 {
402 p.sadd(beta, 1., r);
403 p.add(-beta * omega, v);
404 }
405
406 preconditioner.vmult(y, p);
407 A.vmult(v, y);
408 rhobar = rbar * v;
409 if (std::fabs(rhobar) < additional_data.breakdown)
410 {
411 return IterationResult(true, state, step, res);
412 }
413
414 alpha = rho / rhobar;
415
416 res = std::sqrt(r.add_and_dot(-alpha, v, r));
417
418 // check for early success, see the lac/bicgstab_early testcase as to
419 // why this is necessary
420 //
421 // note: the vector *Vx we pass to the iteration_status signal here is
422 // only the current approximation, not the one we will return with, which
423 // will be x=*Vx + alpha*y
424 if (this->iteration_status(step, res, *Vx) == SolverControl::success)
425 {
426 Vx->add(alpha, y);
427 print_vectors(step, *Vx, r, y);
428 return IterationResult(false, SolverControl::success, step, res);
429 }
430
431 preconditioner.vmult(z, r);
432 A.vmult(t, z);
433 rhobar = t * r;
434 auto t_squared = t * t;
435 if (t_squared < additional_data.breakdown)
436 {
437 return IterationResult(true, state, step, res);
438 }
439 omega = rhobar / (t * t);
440 Vx->add(alpha, y, omega, z);
441
442 if (additional_data.exact_residual)
443 {
444 r.add(-omega, t);
445 res = criterion(A, *Vx, *Vb);
446 }
447 else
448 res = std::sqrt(r.add_and_dot(-omega, t, r));
449
450 state = this->iteration_status(step, res, *Vx);
451 print_vectors(step, *Vx, r, y);
452 }
453 while (state == SolverControl::iterate);
454
455 return IterationResult(false, state, step, res);
456}
457
458
459
460template <typename VectorType>
461template <typename MatrixType, typename PreconditionerType>
462void
463SolverBicgstab<VectorType>::solve(const MatrixType & A,
464 VectorType & x,
465 const VectorType & b,
466 const PreconditionerType &preconditioner)
467{
468 LogStream::Prefix prefix("Bicgstab");
469
470 // Allocate temporary memory.
471 Vr = typename VectorMemory<VectorType>::Pointer(this->memory);
472 Vrbar = typename VectorMemory<VectorType>::Pointer(this->memory);
473 Vp = typename VectorMemory<VectorType>::Pointer(this->memory);
474 Vy = typename VectorMemory<VectorType>::Pointer(this->memory);
475 Vz = typename VectorMemory<VectorType>::Pointer(this->memory);
476 Vt = typename VectorMemory<VectorType>::Pointer(this->memory);
477 Vv = typename VectorMemory<VectorType>::Pointer(this->memory);
478
479 Vr->reinit(x, true);
480 Vrbar->reinit(x, true);
481 Vp->reinit(x, true);
482 Vy->reinit(x, true);
483 Vz->reinit(x, true);
484 Vt->reinit(x, true);
485 Vv->reinit(x, true);
486
487 Vx = &x;
488 Vb = &b;
489
490 step = 0;
491
492 IterationResult state(false, SolverControl::failure, 0, 0);
493 do
494 {
495 state = iterate(A, preconditioner);
496 }
497 while (state.state == SolverControl::iterate);
498
499
500 // Release the temporary memory again.
501 Vr.reset();
502 Vrbar.reset();
503 Vp.reset();
504 Vy.reset();
505 Vz.reset();
506 Vt.reset();
507 Vv.reset();
508
509 // In case of failure: throw exception
510 AssertThrow(state.state == SolverControl::success,
511 SolverControl::NoConvergence(state.last_step,
512 state.last_residual));
513 // Otherwise exit as normal
514}
515
516#endif // DOXYGEN
517
519
520#endif
const VectorType * Vb
VectorMemory< VectorType >::Pointer Vr
VectorMemory< VectorType >::Pointer Vp
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
virtual ~SolverBicgstab() override=default
AdditionalData additional_data
VectorMemory< VectorType >::Pointer Vv
VectorMemory< VectorType >::Pointer Vy
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)
VectorMemory< VectorType >::Pointer Vrbar
VectorMemory< VectorType >::Pointer Vt
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
VectorType * Vx
IterationResult iterate(const MatrixType &A, const PreconditionerType &preconditioner)
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
VectorMemory< VectorType >::Pointer Vz
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
@ failure
Stop iteration, goal not reached.
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
Expression fabs(const Expression &x)
static const char A
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
AdditionalData(const bool exact_residual=true, const double breakdown=std::numeric_limits< typename VectorType::value_type >::min())
IterationResult(const bool breakdown, const SolverControl::State state, const unsigned int last_step, const double last_residual)