Reference documentation for deal.II version 8.5.1
solver_richardson.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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_richardson_h
17 #define dealii__solver_richardson_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 <deal.II/base/subscriptor.h>
25 #include <deal.II/base/signaling_nan.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
31 
61 template <class VectorType = Vector<double> >
62 class SolverRichardson : public Solver<VectorType>
63 {
64 public:
69  {
73  explicit
74  AdditionalData (const double omega = 1,
75  const bool use_preconditioned_residual = false);
76 
80  double omega;
81 
86 
87  };
88 
94  const AdditionalData &data=AdditionalData());
95 
101  const AdditionalData &data=AdditionalData());
102 
106  virtual ~SolverRichardson ();
107 
111  template<typename MatrixType, typename PreconditionerType>
112  void
113  solve (const MatrixType &A,
114  VectorType &x,
115  const VectorType &b,
116  const PreconditionerType &precondition);
117 
121  template<typename MatrixType, typename PreconditionerType>
122  void
123  Tsolve (const MatrixType &A,
124  VectorType &x,
125  const VectorType &b,
126  const PreconditionerType &precondition);
127 
131  void set_omega (const double om=1.);
132 
138  virtual void print_vectors (const unsigned int step,
139  const VectorType &x,
140  const VectorType &r,
141  const VectorType &d) const;
142 
143 protected:
147  virtual typename VectorType::value_type criterion();
148 
153  VectorType *Vr;
159  VectorType *Vd;
160 
165 
172  typename VectorType::value_type res;
173 };
174 
176 /*----------------- Implementation of the Richardson Method ------------------*/
177 
178 #ifndef DOXYGEN
179 
180 template <class VectorType>
181 inline
183 AdditionalData (const double omega,
184  const bool use_preconditioned_residual)
185  :
186  omega(omega),
187  use_preconditioned_residual(use_preconditioned_residual)
188 {}
189 
190 
191 template <class VectorType>
194  const AdditionalData &data)
195  :
196  Solver<VectorType> (cn,mem),
197  Vr(NULL),
198  Vd(NULL),
199  additional_data(data)
200 {}
201 
202 
203 
204 template <class VectorType>
206  const AdditionalData &data)
207  :
208  Solver<VectorType> (cn),
209  Vr(NULL),
210  Vd(NULL),
211  additional_data(data),
212  res(numbers::signaling_nan<double>())
213 {}
214 
215 
216 
217 template <class VectorType>
219 {}
220 
221 
222 template <class VectorType>
223 template <typename MatrixType, typename PreconditionerType>
224 void
225 SolverRichardson<VectorType>::solve (const MatrixType &A,
226  VectorType &x,
227  const VectorType &b,
228  const PreconditionerType &precondition)
229 {
231 
232  double last_criterion = -std::numeric_limits<double>::max();
233 
234  unsigned int iter = 0;
235 
236  // Memory allocation
237  Vr = this->memory.alloc();
238  VectorType &r = *Vr;
239  r.reinit(x);
240  Vd = this->memory.alloc();
241  VectorType &d = *Vd;
242  d.reinit(x);
243 
244  deallog.push("Richardson");
245 
246  try
247  {
248  // Main loop
249  while (conv==SolverControl::iterate)
250  {
251  // Do not use residual,
252  // but do it in 2 steps
253  A.vmult(r,x);
254  r.sadd(-1.,1.,b);
255  precondition.vmult(d,r);
256 
257  // The required norm of the
258  // (preconditioned)
259  // residual is computed in
260  // criterion() and stored
261  // in res.
262  last_criterion = criterion();
263  conv = this->iteration_status (iter, last_criterion, x);
264  if (conv != SolverControl::iterate)
265  break;
266 
267  x.add(additional_data.omega,d);
268  print_vectors(iter,x,r,d);
269 
270  ++iter;
271  }
272  }
273  catch (...)
274  {
275  this->memory.free(Vr);
276  this->memory.free(Vd);
277  deallog.pop();
278  throw;
279  }
280  // Deallocate Memory
281  this->memory.free(Vr);
282  this->memory.free(Vd);
283  deallog.pop();
284 
285  // in case of failure: throw exception
286  if (conv != SolverControl::success)
288  last_criterion));
289  // otherwise exit as normal
290 }
291 
292 
293 template <class VectorType>
294 template <typename MatrixType, typename PreconditionerType>
295 void
296 SolverRichardson<VectorType>::Tsolve (const MatrixType &A,
297  VectorType &x,
298  const VectorType &b,
299  const PreconditionerType &precondition)
300 {
302  double last_criterion = -std::numeric_limits<double>::max();
303 
304  unsigned int iter = 0;
305 
306  // Memory allocation
307  Vr = this->memory.alloc();
308  VectorType &r = *Vr;
309  r.reinit(x);
310  Vd =this-> memory.alloc();
311  VectorType &d = *Vd;
312  d.reinit(x);
313 
314  deallog.push("RichardsonT");
315 
316  try
317  {
318  // Main loop
319  while (conv==SolverControl::iterate)
320  {
321  // Do not use Tresidual,
322  // but do it in 2 steps
323  A.Tvmult(r,x);
324  r.sadd(-1.,1.,b);
325  precondition.Tvmult(d,r);
326 
327  last_criterion = criterion();
328  conv = this->iteration_status (iter, last_criterion, x);
329  if (conv != SolverControl::iterate)
330  break;
331 
332  x.add(additional_data.omega,d);
333  print_vectors(iter,x,r,d);
334 
335  ++iter;
336  }
337  }
338  catch (...)
339  {
340  this->memory.free(Vr);
341  this->memory.free(Vd);
342  deallog.pop();
343  throw;
344  }
345 
346  // Deallocate Memory
347  this->memory.free(Vr);
348  this->memory.free(Vd);
349  deallog.pop();
350  // in case of failure: throw exception
351  if (conv != SolverControl::success)
352  AssertThrow(false, SolverControl::NoConvergence (iter, last_criterion));
353 
354  // otherwise exit as normal
355 }
356 
357 
358 template <class VectorType>
359 void
361  const VectorType &,
362  const VectorType &,
363  const VectorType &) const
364 {}
365 
366 
367 
368 template <class VectorType>
369 inline typename VectorType::value_type
371 {
372  if (!additional_data.use_preconditioned_residual)
373  res = Vr->l2_norm();
374  else
375  res = Vd->l2_norm();
376  return res;
377 }
378 
379 
380 template <class VectorType>
381 inline void
383 {
384  additional_data.omega=om;
385 }
386 
387 #endif // DOXYGEN
388 
389 DEAL_II_NAMESPACE_CLOSE
390 
391 #endif
void pop()
Definition: logstream.cc:306
Continue iteration.
void set_omega(const double om=1.)
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
AdditionalData additional_data
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
SolverRichardson(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Stop iteration, goal reached.
VectorType::value_type res
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
Definition: solver.h:325
void push(const std::string &text)
Definition: logstream.cc:294
void Tsolve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
T signaling_nan()
virtual ~SolverRichardson()
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)
virtual VectorType::value_type criterion()