Reference documentation for deal.II version 9.0.0
solver_bicgstab.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/logstream.h>
22 #include <deal.II/lac/solver.h>
23 #include <deal.II/lac/solver_control.h>
24 #include <cmath>
25 #include <deal.II/base/subscriptor.h>
26 #include <deal.II/base/signaling_nan.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
32 
33 namespace internal
34 {
40  {
41  protected:
45  double alpha;
49  double beta;
53  double omega;
57  double rho;
61  double rhobar;
62 
66  unsigned int step;
67 
71  double res;
72 
78  };
79 }
80 
120 template <typename VectorType = Vector<double> >
121 class SolverBicgstab : public Solver<VectorType>,
123 {
124 public:
136  {
143  explicit
144  AdditionalData(const bool exact_residual = true,
145  const double breakdown = 1.e-10)
148  {}
156  double breakdown;
157  };
158 
164  const AdditionalData &data=AdditionalData());
165 
171  const AdditionalData &data=AdditionalData());
172 
176  virtual ~SolverBicgstab ();
177 
181  template <typename MatrixType, typename PreconditionerType>
182  void
183  solve (const MatrixType &A,
184  VectorType &x,
185  const VectorType &b,
186  const PreconditionerType &preconditioner);
187 
188 protected:
192  VectorType *Vx;
193 
198 
203 
208 
213 
218 
223 
228 
232  const VectorType *Vb;
233 
237  template <typename MatrixType>
238  double criterion (const MatrixType &A, const VectorType &x, const VectorType &b);
239 
245  virtual void print_vectors(const unsigned int step,
246  const VectorType &x,
247  const VectorType &r,
248  const VectorType &d) const;
249 
254 
255 private:
259  template <typename MatrixType>
260  SolverControl::State start(const MatrixType &A);
261 
267  {
268  bool breakdown;
269  SolverControl::State state;
270  unsigned int last_step;
271  double last_residual;
272 
273  IterationResult (const bool breakdown,
274  const SolverControl::State state,
275  const unsigned int last_step,
276  const double last_residual);
277  };
278 
283  template <typename MatrixType, typename PreconditionerType>
285  iterate(const MatrixType &A,
286  const PreconditionerType &preconditioner);
287 };
288 
290 /*-------------------------Inline functions -------------------------------*/
291 
292 #ifndef DOXYGEN
293 
294 
295 template <typename VectorType>
297 (const bool breakdown,
298  const SolverControl::State state,
299  const unsigned int last_step,
300  const double last_residual)
301  :
302  breakdown (breakdown),
303  state (state),
304  last_step (last_step),
305  last_residual (last_residual)
306 {}
307 
308 
309 
310 template <typename VectorType>
313  const AdditionalData &data)
314  :
315  Solver<VectorType>(cn,mem),
316  Vx (nullptr),
317  Vb (nullptr),
318  additional_data(data)
319 {}
320 
321 
322 
323 template <typename VectorType>
325  const AdditionalData &data)
326  :
327  Solver<VectorType>(cn),
328  Vx (nullptr),
329  Vb (nullptr),
330  additional_data(data)
331 {}
332 
333 
334 
335 template <typename VectorType>
337 {}
338 
339 
340 
341 template <typename VectorType>
342 template <typename MatrixType>
343 double
344 SolverBicgstab<VectorType>::criterion (const MatrixType &A,
345  const VectorType &x,
346  const VectorType &b)
347 {
348  A.vmult(*Vt, x);
349  Vt->add(-1.,b);
350  res = Vt->l2_norm();
351 
352  return res;
353 }
354 
355 
356 
357 template <typename VectorType >
358 template <typename MatrixType>
360 SolverBicgstab<VectorType>::start (const MatrixType &A)
361 {
362  A.vmult(*Vr, *Vx);
363  Vr->sadd(-1.,1.,*Vb);
364  res = Vr->l2_norm();
365 
366  return this->iteration_status(step, res, *Vx);
367 }
368 
369 
370 
371 template <typename VectorType>
372 void
374  const VectorType &,
375  const VectorType &,
376  const VectorType &) const
377 {}
378 
379 
380 
381 template <typename VectorType>
382 template <typename MatrixType, typename PreconditionerType>
384 SolverBicgstab<VectorType>::iterate(const MatrixType &A,
385  const PreconditionerType &preconditioner)
386 {
387 //TODO:[GK] Implement "use the length of the computed orthogonal residual" in the BiCGStab method.
389  alpha = omega = rho = 1.;
390 
391  VectorType &r = *Vr;
392  VectorType &rbar = *Vrbar;
393  VectorType &p = *Vp;
394  VectorType &y = *Vy;
395  VectorType &z = *Vz;
396  VectorType &t = *Vt;
397  VectorType &v = *Vv;
398 
399  rbar = r;
400  bool startup = true;
401 
402  do
403  {
404  ++step;
405 
406  rhobar = r*rbar;
407  beta = rhobar * alpha / (rho * omega);
408  rho = rhobar;
409  if (startup == true)
410  {
411  p = r;
412  startup = false;
413  }
414  else
415  {
416  p.sadd(beta, 1., r);
417  p.add(-beta*omega, v);
418  }
419 
420  preconditioner.vmult(y,p);
421  A.vmult(v,y);
422  rhobar = rbar * v;
423 
424  alpha = rho/rhobar;
425 
426 //TODO:[?] Find better breakdown criterion
427 
428  if (std::fabs(alpha) > 1.e10)
429  return IterationResult(true, state, step, res);
430 
431  res = std::sqrt(r.add_and_dot(-alpha, v, r));
432 
433  // check for early success, see the lac/bicgstab_early testcase as to
434  // why this is necessary
435  //
436  // note: the vector *Vx we pass to the iteration_status signal here is only
437  // the current approximation, not the one we will return with,
438  // which will be x=*Vx + alpha*y
440  {
441  Vx->add(alpha, y);
442  print_vectors(step, *Vx, r, y);
443  return IterationResult(false, SolverControl::success, step, res);
444  }
445 
446  preconditioner.vmult(z,r);
447  A.vmult(t,z);
448  rhobar = t*r;
449  omega = rhobar/(t*t);
450  Vx->add(alpha, y, omega, z);
451 
453  {
454  r.add(-omega, t);
455  res = criterion(A, *Vx, *Vb);
456  }
457  else
458  res = std::sqrt(r.add_and_dot(-omega, t, r));
459 
460  state = this->iteration_status(step, res, *Vx);
461  print_vectors(step, *Vx, r, y);
462  }
463  while (state == SolverControl::iterate);
464  return IterationResult(false, state, step, res);
465 }
466 
467 
468 template <typename VectorType>
469 template <typename MatrixType, typename PreconditionerType>
470 void
471 SolverBicgstab<VectorType>::solve(const MatrixType &A,
472  VectorType &x,
473  const VectorType &b,
474  const PreconditionerType &preconditioner)
475 {
476  LogStream::Prefix prefix("Bicgstab");
477  Vr = typename VectorMemory<VectorType>::Pointer(this->memory);
479  Vp = typename VectorMemory<VectorType>::Pointer(this->memory);
480  Vy = typename VectorMemory<VectorType>::Pointer(this->memory);
481  Vz = typename VectorMemory<VectorType>::Pointer(this->memory);
482  Vt = typename VectorMemory<VectorType>::Pointer(this->memory);
483  Vv = typename VectorMemory<VectorType>::Pointer(this->memory);
484 
485  Vr->reinit(x, true);
486  Vrbar->reinit(x, true);
487  Vp->reinit(x, true);
488  Vy->reinit(x, true);
489  Vz->reinit(x, true);
490  Vt->reinit(x, true);
491  Vv->reinit(x, true);
492 
493  Vx = &x;
494  Vb = &b;
495 
496  step = 0;
497 
498  IterationResult state(false,SolverControl::failure,0,0);
499 
500  // iterate while the inner iteration returns a breakdown
501  do
502  {
503  if (step != 0)
504  deallog << "Restart step " << step << std::endl;
505  if (start(A) == SolverControl::success)
506  {
507  state.state = SolverControl::success;
508  break;
509  }
510  state = iterate(A, preconditioner);
511  ++step;
512  }
513  while (state.breakdown == true);
514 
515  // in case of failure: throw exception
516  AssertThrow(state.state == SolverControl::success,
517  SolverControl::NoConvergence (state.last_step,
518  state.last_residual));
519  // otherwise exit as normal
520 }
521 
522 #endif // DOXYGEN
523 
524 DEAL_II_NAMESPACE_CLOSE
525 
526 #endif
Stop iteration, goal not reached.
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Continue iteration.
VectorMemory< VectorType >::Pointer Vrbar
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
virtual ~SolverBicgstab()
VectorMemory< VectorType >::Pointer Vv
VectorMemory< VectorType >::Pointer Vp
VectorMemory< VectorType >::Pointer Vr
Stop iteration, goal reached.
AdditionalData additional_data
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
VectorMemory< VectorType >::Pointer Vt
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorMemory< VectorType > & memory
Definition: solver.h:399
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate), StateCombiner > iteration_status
Definition: solver.h:447
AdditionalData(const bool exact_residual=true, const double breakdown=1.e-10)
Definition: solver.h:322
const VectorType * Vb
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
IterationResult iterate(const MatrixType &A, const PreconditionerType &preconditioner)
VectorMemory< VectorType >::Pointer Vy
VectorMemory< VectorType >::Pointer Vz
SolverControl::State start(const MatrixType &A)
VectorType * Vx
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)