Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
solver_minres.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 2022 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.md at
12// the top level directory of deal.II.
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
25
26#include <deal.II/lac/solver.h>
28
29#include <cmath>
30
32
35
67template <class VectorType = Vector<double>>
68class SolverMinRes : public SolverBase<VectorType>
69{
70public:
76 {};
77
83 const AdditionalData & data = AdditionalData());
84
90 const AdditionalData &data = AdditionalData());
91
95 virtual ~SolverMinRes() override = default;
96
100 template <typename MatrixType, typename PreconditionerType>
101 void
102 solve(const MatrixType & A,
103 VectorType & x,
104 const VectorType & b,
105 const PreconditionerType &preconditioner);
106
117
118protected:
122 virtual double
124
130 virtual void
131 print_vectors(const unsigned int step,
132 const VectorType & x,
133 const VectorType & r,
134 const VectorType & d) const;
135
142 double res2;
143};
144
146/*------------------------- Implementation ----------------------------*/
147
148#ifndef DOXYGEN
149
150template <class VectorType>
153 const AdditionalData &)
154 : SolverBase<VectorType>(cn, mem)
155 , res2(numbers::signaling_nan<double>())
156{}
157
158
159
160template <class VectorType>
162 const AdditionalData &)
163 : SolverBase<VectorType>(cn)
164 , res2(numbers::signaling_nan<double>())
165{}
166
167
168
169template <class VectorType>
170double
172{
173 return res2;
174}
175
176
177template <class VectorType>
178void
180 const VectorType &,
181 const VectorType &,
182 const VectorType &) const
183{}
184
185
186
187template <class VectorType>
188template <typename MatrixType, typename PreconditionerType>
189void
190SolverMinRes<VectorType>::solve(const MatrixType & A,
191 VectorType & x,
192 const VectorType & b,
193 const PreconditionerType &preconditioner)
194{
195 LogStream::Prefix prefix("minres");
196
197 // Memory allocation
198 typename VectorMemory<VectorType>::Pointer Vu0(this->memory);
199 typename VectorMemory<VectorType>::Pointer Vu1(this->memory);
200 typename VectorMemory<VectorType>::Pointer Vu2(this->memory);
201
202 typename VectorMemory<VectorType>::Pointer Vm0(this->memory);
203 typename VectorMemory<VectorType>::Pointer Vm1(this->memory);
204 typename VectorMemory<VectorType>::Pointer Vm2(this->memory);
205
206 typename VectorMemory<VectorType>::Pointer Vv(this->memory);
207
208 // define some aliases for simpler access
209 using vecptr = VectorType *;
210 vecptr u[3] = {Vu0.get(), Vu1.get(), Vu2.get()};
211 vecptr m[3] = {Vm0.get(), Vm1.get(), Vm2.get()};
212 VectorType &v = *Vv;
213
214 // resize the vectors, but do not set the values since they'd be overwritten
215 // soon anyway.
216 u[0]->reinit(b, true);
217 u[1]->reinit(b, true);
218 u[2]->reinit(b, true);
219 m[0]->reinit(b, true);
220 m[1]->reinit(b, true);
221 m[2]->reinit(b, true);
222 v.reinit(b, true);
223
224 // some values needed
225 double delta[3] = {0, 0, 0};
226 double f[2] = {0, 0};
227 double e[2] = {0, 0};
228
229 double r_l2 = 0;
230 double r0 = 0;
231 double tau = 0;
232 double c = 0;
233 double s = 0;
234 double d_ = 0;
235
236 // The iteration step.
237 unsigned int j = 1;
238
239
240 // Start of the solution process
241 A.vmult(*m[0], x);
242 *u[1] = b;
243 *u[1] -= *m[0];
244 // Precondition is applied.
245 // The preconditioner has to be
246 // positive definite and symmetric
247
248 // M v = u[1]
249 preconditioner.vmult(v, *u[1]);
250
251 delta[1] = v * (*u[1]);
252 // Preconditioner positive
253 Assert(delta[1] >= 0, ExcPreconditionerNotDefinite());
254
255 r0 = std::sqrt(delta[1]);
256 r_l2 = r0;
257
258
259 u[0]->reinit(b);
260 delta[0] = 1.;
261 m[0]->reinit(b);
262 m[1]->reinit(b);
263 m[2]->reinit(b);
264
265 SolverControl::State conv = this->iteration_status(0, r_l2, x);
266 while (conv == SolverControl::iterate)
267 {
268 if (delta[1] != 0)
269 v *= 1. / std::sqrt(delta[1]);
270 else
271 v.reinit(b);
272
273 A.vmult(*u[2], v);
274 u[2]->add(-std::sqrt(delta[1] / delta[0]), *u[0]);
275
276 const double gamma = *u[2] * v;
277 u[2]->add(-gamma / std::sqrt(delta[1]), *u[1]);
278 *m[0] = v;
279
280 // precondition: solve M v = u[2]
281 // Preconditioner has to be positive
282 // definite and symmetric.
283 preconditioner.vmult(v, *u[2]);
284
285 delta[2] = v * (*u[2]);
286
287 Assert(delta[2] >= 0, ExcPreconditionerNotDefinite());
288
289 if (j == 1)
290 {
291 d_ = gamma;
292 e[1] = std::sqrt(delta[2]);
293 }
294 if (j > 1)
295 {
296 d_ = s * e[0] - c * gamma;
297 e[0] = c * e[0] + s * gamma;
298 f[1] = s * std::sqrt(delta[2]);
299 e[1] = -c * std::sqrt(delta[2]);
300 }
301
302 const double d = std::sqrt(d_ * d_ + delta[2]);
303
304 if (j > 1)
305 tau *= s / c;
306 c = d_ / d;
307 tau *= c;
308
309 s = std::sqrt(delta[2]) / d;
310
311 if (j == 1)
312 tau = r0 * c;
313
314 m[0]->add(-e[0], *m[1]);
315 if (j > 1)
316 m[0]->add(-f[0], *m[2]);
317 *m[0] *= 1. / d;
318 x.add(tau, *m[0]);
319 r_l2 *= std::fabs(s);
320
321 conv = this->iteration_status(j, r_l2, x);
322
323 // next iteration step
324 ++j;
325 // All vectors have to be shifted
326 // one iteration step.
327 // This should be changed one time.
328 swap(*m[2], *m[1]);
329 swap(*m[1], *m[0]);
330
331 // likewise, but reverse direction:
332 // u[0] = u[1];
333 // u[1] = u[2];
334 swap(*u[0], *u[1]);
335 swap(*u[1], *u[2]);
336
337 // these are scalars, so need
338 // to bother
339 f[0] = f[1];
340 e[0] = e[1];
341 delta[0] = delta[1];
342 delta[1] = delta[2];
343 }
344
345 // in case of failure: throw exception
348
349 // otherwise exit as normal
350}
351
352#endif // DOXYGEN
353
355
356#endif
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
virtual ~SolverMinRes() override=default
SolverMinRes(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverMinRes(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcPreconditionerNotDefinite()
#define Assert(cond, exc)
Definition: exceptions.h:1473
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
virtual double criterion()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
long double gamma(const unsigned int n)
T signaling_nan()
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
void swap(SmartPointer< T, P > &t1, SmartPointer< T, Q > &t2)
Definition: smartpointer.h:447