Reference documentation for deal.II version 8.5.1
solver_bicgstab.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 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 
72 template <typename VectorType = Vector<double> >
73 class SolverBicgstab : public Solver<VectorType>
74 {
75 public:
87  {
94  explicit
95  AdditionalData(const bool exact_residual = true,
96  const double breakdown = 1.e-10)
99  {}
107  double breakdown;
108  };
109 
115  const AdditionalData &data=AdditionalData());
116 
122  const AdditionalData &data=AdditionalData());
123 
127  virtual ~SolverBicgstab ();
128 
132  template<typename MatrixType, typename PreconditionerType>
133  void
134  solve (const MatrixType &A,
135  VectorType &x,
136  const VectorType &b,
137  const PreconditionerType &precondition);
138 
139 protected:
143  template <typename MatrixType>
144  double criterion (const MatrixType &A, const VectorType &x, const VectorType &b);
145 
151  virtual void print_vectors(const unsigned int step,
152  const VectorType &x,
153  const VectorType &r,
154  const VectorType &d) const;
155 
159  VectorType *Vx;
163  VectorType *Vr;
167  VectorType *Vrbar;
171  VectorType *Vp;
175  VectorType *Vy;
179  VectorType *Vz;
183  VectorType *Vt;
187  VectorType *Vv;
191  const VectorType *Vb;
192 
196  double alpha;
200  double beta;
204  double omega;
208  double rho;
212  double rhobar;
213 
217  unsigned int step;
218 
222  double res;
223 
228 
229 private:
233  template <typename MatrixType>
234  SolverControl::State start(const MatrixType &A);
235 
241  {
242  bool breakdown;
243  SolverControl::State state;
244  unsigned int last_step;
245  double last_residual;
246 
247  IterationResult (const bool breakdown,
248  const SolverControl::State state,
249  const unsigned int last_step,
250  const double last_residual);
251  };
252 
257  template<typename MatrixType, typename PreconditionerType>
259  iterate(const MatrixType &A,
260  const PreconditionerType &precondition);
261 };
262 
264 /*-------------------------Inline functions -------------------------------*/
265 
266 #ifndef DOXYGEN
267 
268 
269 template<typename VectorType>
271 (const bool breakdown,
272  const SolverControl::State state,
273  const unsigned int last_step,
274  const double last_residual)
275  :
276  breakdown (breakdown),
277  state (state),
278  last_step (last_step),
279  last_residual (last_residual)
280 {}
281 
282 
283 template<typename VectorType>
286  const AdditionalData &data)
287  :
288  Solver<VectorType>(cn,mem),
289  additional_data(data)
290 {}
291 
292 
293 
294 template<typename VectorType>
296  const AdditionalData &data)
297  :
298  Solver<VectorType>(cn),
299  Vx(NULL),
300  Vr(NULL),
301  Vrbar(NULL),
302  Vp(NULL),
303  Vy(NULL),
304  Vz(NULL),
305  Vt(NULL),
306  Vv(NULL),
307  Vb(NULL),
308  alpha(0.),
309  beta(0.),
310  omega(0.),
311  rho(0.),
312  rhobar(0.),
313  step(numbers::invalid_unsigned_int),
314  res(numbers::signaling_nan<double>()),
315  additional_data(data)
316 {}
317 
318 
319 
320 template<typename VectorType>
322 {}
323 
324 
325 
326 template <typename VectorType>
327 template <typename MatrixType>
328 double
329 SolverBicgstab<VectorType>::criterion (const MatrixType &A,
330  const VectorType &x,
331  const VectorType &b)
332 {
333  A.vmult(*Vt, x);
334  Vt->add(-1.,b);
335  res = Vt->l2_norm();
336 
337  return res;
338 }
339 
340 
341 
342 template <typename VectorType >
343 template <typename MatrixType>
345 SolverBicgstab<VectorType>::start (const MatrixType &A)
346 {
347  A.vmult(*Vr, *Vx);
348  Vr->sadd(-1.,1.,*Vb);
349  res = Vr->l2_norm();
350 
351  return this->iteration_status(step, res, *Vx);
352 }
353 
354 
355 
356 template<typename VectorType>
357 void
359  const VectorType &,
360  const VectorType &,
361  const VectorType &) const
362 {}
363 
364 
365 
366 template<typename VectorType>
367 template<typename MatrixType, typename PreconditionerType>
369 SolverBicgstab<VectorType>::iterate(const MatrixType &A,
370  const PreconditionerType &precondition)
371 {
372 //TODO:[GK] Implement "use the length of the computed orthogonal residual" in the BiCGStab method.
374  alpha = omega = rho = 1.;
375 
376  VectorType &r = *Vr;
377  VectorType &rbar = *Vrbar;
378  VectorType &p = *Vp;
379  VectorType &y = *Vy;
380  VectorType &z = *Vz;
381  VectorType &t = *Vt;
382  VectorType &v = *Vv;
383 
384  rbar = r;
385  bool startup = true;
386 
387  do
388  {
389  ++step;
390 
391  rhobar = r*rbar;
392  beta = rhobar * alpha / (rho * omega);
393  rho = rhobar;
394  if (startup == true)
395  {
396  p = r;
397  startup = false;
398  }
399  else
400  {
401  p.sadd(beta, 1., r);
402  p.add(-beta*omega, v);
403  }
404 
405  precondition.vmult(y,p);
406  A.vmult(v,y);
407  rhobar = rbar * v;
408 
409  alpha = rho/rhobar;
410 
411 //TODO:[?] Find better breakdown criterion
412 
413  if (std::fabs(alpha) > 1.e10)
414  return IterationResult(true, state, step, res);
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 only
422  // the current approximation, not the one we will return with,
423  // which will be x=*Vx + alpha*y
425  {
426  Vx->add(alpha, y);
427  print_vectors(step, *Vx, r, y);
428  return IterationResult(false, SolverControl::success, step, res);
429  }
430 
431  precondition.vmult(z,r);
432  A.vmult(t,z);
433  rhobar = t*r;
434  omega = rhobar/(t*t);
435  Vx->add(alpha, y, omega, z);
436 
438  {
439  r.add(-omega, t);
440  res = criterion(A, *Vx, *Vb);
441  }
442  else
443  res = std::sqrt(r.add_and_dot(-omega, t, r));
444 
445  state = this->iteration_status(step, res, *Vx);
446  print_vectors(step, *Vx, r, y);
447  }
448  while (state == SolverControl::iterate);
449  return IterationResult(false, state, step, res);
450 }
451 
452 
453 template<typename VectorType>
454 template<typename MatrixType, typename PreconditionerType>
455 void
456 SolverBicgstab<VectorType>::solve(const MatrixType &A,
457  VectorType &x,
458  const VectorType &b,
459  const PreconditionerType &precondition)
460 {
461  deallog.push("Bicgstab");
462  Vr = this->memory.alloc();
463  Vr->reinit(x, true);
464  Vrbar = this->memory.alloc();
465  Vrbar->reinit(x, true);
466  Vp = this->memory.alloc();
467  Vp->reinit(x, true);
468  Vy = this->memory.alloc();
469  Vy->reinit(x, true);
470  Vz = this->memory.alloc();
471  Vz->reinit(x, true);
472  Vt = this->memory.alloc();
473  Vt->reinit(x, true);
474  Vv = this->memory.alloc();
475  Vv->reinit(x, true);
476 
477  Vx = &x;
478  Vb = &b;
479 
480  step = 0;
481 
482  IterationResult state(false,SolverControl::failure,0,0);
483 
484  // iterate while the inner iteration returns a breakdown
485  do
486  {
487  if (step != 0)
488  deallog << "Restart step " << step << std::endl;
489  if (start(A) == SolverControl::success)
490  {
491  state.state = SolverControl::success;
492  break;
493  }
494  state = iterate(A, precondition);
495  ++step;
496  }
497  while (state.breakdown == true);
498 
499  this->memory.free(Vr);
500  this->memory.free(Vrbar);
501  this->memory.free(Vp);
502  this->memory.free(Vy);
503  this->memory.free(Vz);
504  this->memory.free(Vt);
505  this->memory.free(Vv);
506 
507  deallog.pop();
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 
518 DEAL_II_NAMESPACE_CLOSE
519 
520 #endif
Stop iteration, goal not reached.
SolverBicgstab(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
VectorType * Vp
void pop()
Definition: logstream.cc:306
IterationResult iterate(const MatrixType &A, const PreconditionerType &precondition)
Continue iteration.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual ~SolverBicgstab()
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
unsigned int step
Stop iteration, goal reached.
AdditionalData additional_data
VectorType * Vv
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType * Vz
VectorMemory< VectorType > & memory
Definition: solver.h:402
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate), StateCombiner > iteration_status
Definition: solver.h:450
AdditionalData(const bool exact_residual=true, const double breakdown=1.e-10)
Definition: solver.h:325
const VectorType * Vb
void push(const std::string &text)
Definition: logstream.cc:294
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
T signaling_nan()
VectorType * Vr
SolverControl::State start(const MatrixType &A)
VectorType * Vx
VectorType * Vt
VectorType * Vrbar
VectorType * Vy
double criterion(const MatrixType &A, const VectorType &x, const VectorType &b)