Reference documentation for deal.II version 8.5.1
solver_minres.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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_minres_h
17 #define dealii__solver_minres_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/lac/solver.h>
22 #include <deal.II/lac/solver_control.h>
23 #include <deal.II/base/logstream.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 
68 template <class VectorType = Vector<double> >
69 class SolverMinRes : public Solver<VectorType>
70 {
71 public:
77  {
78  };
79 
85  const AdditionalData &data=AdditionalData());
86 
92  const AdditionalData &data=AdditionalData());
93 
97  virtual ~SolverMinRes ();
98 
102  template<typename MatrixType, typename PreconditionerType>
103  void
104  solve (const MatrixType &A,
105  VectorType &x,
106  const VectorType &b,
107  const PreconditionerType &precondition);
108 
119 
120 protected:
124  virtual double criterion();
130  virtual void print_vectors(const unsigned int step,
131  const VectorType &x,
132  const VectorType &r,
133  const VectorType &d) const;
134 
139  VectorType *Vu0, *Vu1, *Vu2;
140  VectorType *Vm0, *Vm1, *Vm2;
141  VectorType *Vv;
142 
149  double res2;
150 };
151 
153 /*------------------------- Implementation ----------------------------*/
154 
155 #ifndef DOXYGEN
156 
157 template<class VectorType>
160  const AdditionalData &)
161  :
162  Solver<VectorType>(cn,mem)
163 {}
164 
165 
166 
167 template<class VectorType>
169  const AdditionalData &)
170  :
171  Solver<VectorType>(cn),
172  Vu0(NULL),
173  Vu1(NULL),
174  Vu2(NULL),
175  Vm0(NULL),
176  Vm1(NULL),
177  Vm2(NULL),
178  Vv(NULL),
179  res2(numbers::signaling_nan<double>())
180 {}
181 
182 
183 template<class VectorType>
185 {}
186 
187 
188 
189 template<class VectorType>
190 double
192 {
193  return res2;
194 }
195 
196 
197 template<class VectorType>
198 void
199 SolverMinRes<VectorType>::print_vectors(const unsigned int,
200  const VectorType &,
201  const VectorType &,
202  const VectorType &) const
203 {}
204 
205 
206 
207 template<class VectorType>
208 template<typename MatrixType, typename PreconditionerType>
209 void
210 SolverMinRes<VectorType>::solve (const MatrixType &A,
211  VectorType &x,
212  const VectorType &b,
213  const PreconditionerType &precondition)
214 {
215  deallog.push("minres");
216 
217  // Memory allocation
218  Vu0 = this->memory.alloc();
219  Vu1 = this->memory.alloc();
220  Vu2 = this->memory.alloc();
221  Vv = this->memory.alloc();
222  Vm0 = this->memory.alloc();
223  Vm1 = this->memory.alloc();
224  Vm2 = this->memory.alloc();
225  // define some aliases for simpler access
226  typedef VectorType *vecptr;
227  vecptr u[3] = {Vu0, Vu1, Vu2};
228  vecptr m[3] = {Vm0, Vm1, Vm2};
229  VectorType &v = *Vv;
230  // resize the vectors, but do not set
231  // the values since they'd be overwritten
232  // soon anyway.
233  u[0]->reinit(b,true);
234  u[1]->reinit(b,true);
235  u[2]->reinit(b,true);
236  m[0]->reinit(b,true);
237  m[1]->reinit(b,true);
238  m[2]->reinit(b,true);
239  v.reinit(b,true);
240 
241  // some values needed
242  double delta[3] = { 0, 0, 0 };
243  double f[2] = { 0, 0 };
244  double e[2] = { 0, 0 };
245 
246  double r_l2 = 0;
247  double r0 = 0;
248  double tau = 0;
249  double c = 0;
250  double s = 0;
251  double d_ = 0;
252 
253  // The iteration step.
254  unsigned int j = 1;
255 
256 
257  // Start of the solution process
258  A.vmult(*m[0],x);
259  *u[1] = b;
260  *u[1] -= *m[0];
261  // Precondition is applied.
262  // The preconditioner has to be
263  // positive definite and symmetric
264 
265  // M v = u[1]
266  precondition.vmult (v,*u[1]);
267 
268  delta[1] = v * (*u[1]);
269  // Preconditioner positive
270  Assert (delta[1]>=0, ExcPreconditionerNotDefinite());
271 
272  r0 = std::sqrt(delta[1]);
273  r_l2 = r0;
274 
275 
276  u[0]->reinit(b);
277  delta[0] = 1.;
278  m[0]->reinit(b);
279  m[1]->reinit(b);
280  m[2]->reinit(b);
281 
282  SolverControl::State conv = this->iteration_status(0, r_l2, x);
283  while (conv==SolverControl::iterate)
284  {
285  if (delta[1]!=0)
286  v *= 1./std::sqrt(delta[1]);
287  else
288  v.reinit(b);
289 
290  A.vmult(*u[2],v);
291  u[2]->add (-std::sqrt(delta[1]/delta[0]), *u[0]);
292 
293  const double gamma = *u[2] * v;
294  u[2]->add (-gamma / std::sqrt(delta[1]), *u[1]);
295  *m[0] = v;
296 
297  // precondition: solve M v = u[2]
298  // Preconditioner has to be positive
299  // definite and symmetric.
300  precondition.vmult(v,*u[2]);
301 
302  delta[2] = v * (*u[2]);
303 
304  Assert (delta[2]>=0, ExcPreconditionerNotDefinite());
305 
306  if (j==1)
307  {
308  d_ = gamma;
309  e[1] = std::sqrt(delta[2]);
310  }
311  if (j>1)
312  {
313  d_ = s * e[0] - c * gamma;
314  e[0] = c * e[0] + s * gamma;
315  f[1] = s * std::sqrt(delta[2]);
316  e[1] = -c * std::sqrt(delta[2]);
317  }
318 
319  const double d = std::sqrt (d_*d_ + delta[2]);
320 
321  if (j>1)
322  tau *= s / c;
323  c = d_ / d;
324  tau *= c;
325 
326  s = std::sqrt(delta[2]) / d;
327 
328  if (j==1)
329  tau = r0 * c;
330 
331  m[0]->add (-e[0], *m[1]);
332  if (j>1)
333  m[0]->add (-f[0], *m[2]);
334  *m[0] *= 1./d;
335  x.add (tau, *m[0]);
336  r_l2 *= std::fabs(s);
337 
338  conv = this->iteration_status(j, r_l2, x);
339 
340  // next iteration step
341  ++j;
342  // All vectors have to be shifted
343  // one iteration step.
344  // This should be changed one time.
345  //
346  // the previous code was like this:
347  // m[2] = m[1];
348  // m[1] = m[0];
349  // but it can be made more efficient,
350  // since the value of m[0] is no more
351  // needed in the next iteration
352  swap (*m[2], *m[1]);
353  swap (*m[1], *m[0]);
354 
355  // likewise, but reverse direction:
356  // u[0] = u[1];
357  // u[1] = u[2];
358  swap (*u[0], *u[1]);
359  swap (*u[1], *u[2]);
360 
361  // these are scalars, so need
362  // to bother
363  f[0] = f[1];
364  e[0] = e[1];
365  delta[0] = delta[1];
366  delta[1] = delta[2];
367  }
368 
369  // Deallocation of Memory
370  this->memory.free(Vu0);
371  this->memory.free(Vu1);
372  this->memory.free(Vu2);
373  this->memory.free(Vv);
374  this->memory.free(Vm0);
375  this->memory.free(Vm1);
376  this->memory.free(Vm2);
377  // Output
378  deallog.pop ();
379 
380  // in case of failure: throw exception
382  SolverControl::NoConvergence (j, r_l2));
383 
384  // otherwise exit as normal
385 }
386 
387 #endif // DOXYGEN
388 
389 DEAL_II_NAMESPACE_CLOSE
390 
391 #endif
void pop()
Definition: logstream.cc:306
Continue iteration.
SolverMinRes(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
VectorType * Vu0
void swap(BlockIndices &u, BlockIndices &v)
virtual ~SolverMinRes()
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:313
#define DeclException0(Exception0)
Definition: exceptions.h:541
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Definition: solver.h:325
void push(const std::string &text)
Definition: logstream.cc:294
T signaling_nan()
virtual double criterion()
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
static ::ExceptionBase & ExcPreconditionerNotDefinite()
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const