Reference documentation for deal.II version 9.0.0
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/signaling_nan.h>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
30 
60 template <class VectorType = Vector<double> >
61 class SolverRichardson : public Solver<VectorType>
62 {
63 public:
68  {
72  explicit
73  AdditionalData (const double omega = 1,
74  const bool use_preconditioned_residual = false);
75 
79  double omega;
80 
85 
86  };
87 
93  const AdditionalData &data=AdditionalData());
94 
100  const AdditionalData &data=AdditionalData());
101 
105  virtual ~SolverRichardson () = default;
106 
110  template <typename MatrixType, typename PreconditionerType>
111  void
112  solve (const MatrixType &A,
113  VectorType &x,
114  const VectorType &b,
115  const PreconditionerType &preconditioner);
116 
120  template <typename MatrixType, typename PreconditionerType>
121  void
122  Tsolve (const MatrixType &A,
123  VectorType &x,
124  const VectorType &b,
125  const PreconditionerType &preconditioner);
126 
130  void set_omega (const double om=1.);
131 
137  virtual void print_vectors (const unsigned int step,
138  const VectorType &x,
139  const VectorType &r,
140  const VectorType &d) const;
141 
142 protected:
149  virtual typename VectorType::value_type criterion(const VectorType &r,
150  const VectorType &d) const;
151 
156 };
157 
159 /*----------------- Implementation of the Richardson Method ------------------*/
160 
161 #ifndef DOXYGEN
162 
163 template <class VectorType>
164 inline
166 AdditionalData (const double omega,
167  const bool use_preconditioned_residual)
168  :
169  omega(omega),
170  use_preconditioned_residual(use_preconditioned_residual)
171 {}
172 
173 
174 template <class VectorType>
177  const AdditionalData &data)
178  :
179  Solver<VectorType> (cn,mem),
180  additional_data(data)
181 {}
182 
183 
184 
185 template <class VectorType>
187  const AdditionalData &data)
188  :
189  Solver<VectorType> (cn),
190  additional_data(data)
191 {}
192 
193 
194 
195 template <class VectorType>
196 template <typename MatrixType, typename PreconditionerType>
197 void
198 SolverRichardson<VectorType>::solve (const MatrixType &A,
199  VectorType &x,
200  const VectorType &b,
201  const PreconditionerType &preconditioner)
202 {
204 
205  double last_criterion = -std::numeric_limits<double>::max();
206 
207  unsigned int iter = 0;
208 
209  // Memory allocation.
210  // 'Vr' holds the residual, 'Vd' the preconditioned residual
211  typename VectorMemory<VectorType>::Pointer Vr (this->memory);
212  typename VectorMemory<VectorType>::Pointer Vd (this->memory);
213 
214  VectorType &r = *Vr;
215  r.reinit(x);
216 
217  VectorType &d = *Vd;
218  d.reinit(x);
219 
220  LogStream::Prefix prefix("Richardson");
221 
222  // Main loop
223  while (conv==SolverControl::iterate)
224  {
225  // Do not use residual,
226  // but do it in 2 steps
227  A.vmult(r,x);
228  r.sadd(-1.,1.,b);
229  preconditioner.vmult(d,r);
230 
231  // get the required norm of the (possibly preconditioned)
232  // residual
233  last_criterion = criterion(r, d);
234  conv = this->iteration_status (iter, last_criterion, x);
235  if (conv != SolverControl::iterate)
236  break;
237 
238  x.add(additional_data.omega,d);
239  print_vectors(iter,x,r,d);
240 
241  ++iter;
242  }
243 
244  // in case of failure: throw exception
245  if (conv != SolverControl::success)
247  last_criterion));
248  // otherwise exit as normal
249 }
250 
251 
252 
253 template <class VectorType>
254 template <typename MatrixType, typename PreconditionerType>
255 void
256 SolverRichardson<VectorType>::Tsolve (const MatrixType &A,
257  VectorType &x,
258  const VectorType &b,
259  const PreconditionerType &preconditioner)
260 {
262  double last_criterion = -std::numeric_limits<double>::max();
263 
264  unsigned int iter = 0;
265 
266  // Memory allocation.
267  // 'Vr' holds the residual, 'Vd' the preconditioned residual
268  typename VectorMemory<VectorType>::Pointer Vr (this->memory);
269  typename VectorMemory<VectorType>::Pointer Vd (this->memory);
270 
271  VectorType &r = *Vr;
272  r.reinit(x);
273 
274  VectorType &d = *Vd;
275  d.reinit(x);
276 
277  LogStream::Prefix prefix("RichardsonT");
278 
279  // Main loop
280  while (conv==SolverControl::iterate)
281  {
282  // Do not use Tresidual,
283  // but do it in 2 steps
284  A.Tvmult(r,x);
285  r.sadd(-1.,1.,b);
286  preconditioner.Tvmult(d,r);
287 
288  last_criterion = criterion(r, d);
289  conv = this->iteration_status (iter, last_criterion, x);
290  if (conv != SolverControl::iterate)
291  break;
292 
293  x.add(additional_data.omega,d);
294  print_vectors(iter,x,r,d);
295 
296  ++iter;
297  }
298 
299  // in case of failure: throw exception
300  if (conv != SolverControl::success)
301  AssertThrow(false, SolverControl::NoConvergence (iter, last_criterion));
302 
303  // otherwise exit as normal
304 }
305 
306 
307 template <class VectorType>
308 void
310  const VectorType &,
311  const VectorType &,
312  const VectorType &) const
313 {}
314 
315 
316 
317 template <class VectorType>
318 inline
319 typename VectorType::value_type
320 SolverRichardson<VectorType>::criterion(const VectorType &r,
321  const VectorType &d) const
322 {
323  if (!additional_data.use_preconditioned_residual)
324  return r.l2_norm();
325  else
326  return d.l2_norm();
327 }
328 
329 
330 template <class VectorType>
331 inline void
333 {
334  additional_data.omega=om;
335 }
336 
337 #endif // DOXYGEN
338 
339 DEAL_II_NAMESPACE_CLOSE
340 
341 #endif
virtual ~SolverRichardson()=default
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
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
void Tsolve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
SolverRichardson(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Stop iteration, goal reached.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Definition: solver.h:322
virtual VectorType::value_type criterion(const VectorType &r, const VectorType &d) const
AdditionalData(const double omega=1, const bool use_preconditioned_residual=false)