Reference documentation for deal.II version 9.3.3
\(\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\}}\)
solver_idr.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 2020 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_idr_h
17#define dealii_solver_idr_h
18
19
20#include <deal.II/base/config.h>
21
26
28#include <deal.II/lac/solver.h>
30
31#include <cmath>
32#include <random>
33
35
38
39namespace internal
40{
44 namespace SolverIDRImplementation
45 {
50 template <typename VectorType>
52 {
53 public:
57 TmpVectors(const unsigned int s_param, VectorMemory<VectorType> &vmem);
58
62 ~TmpVectors() = default;
63
68 VectorType &operator[](const unsigned int i) const;
69
76 VectorType &
77 operator()(const unsigned int i, const VectorType &temp);
78
79 private:
84
88 std::vector<typename VectorMemory<VectorType>::Pointer> data;
89 };
90 } // namespace SolverIDRImplementation
91} // namespace internal
92
115template <class VectorType = Vector<double>>
116class SolverIDR : public SolverBase<VectorType>
117{
118public:
123 {
127 explicit AdditionalData(const unsigned int s = 2)
128 : s(s)
129 {}
130
131 const unsigned int s;
132 };
133
139 const AdditionalData & data = AdditionalData());
140
146 const AdditionalData &data = AdditionalData());
147
151 virtual ~SolverIDR() override = default;
152
156 template <typename MatrixType, typename PreconditionerType>
157 void
158 solve(const MatrixType & A,
159 VectorType & x,
160 const VectorType & b,
161 const PreconditionerType &preconditioner);
162
163protected:
169 virtual void
170 print_vectors(const unsigned int step,
171 const VectorType & x,
172 const VectorType & r,
173 const VectorType & d) const;
174
175private:
180};
181
183/*------------------------- Implementation ----------------------------*/
184
185#ifndef DOXYGEN
186
187
188namespace internal
189{
190 namespace SolverIDRImplementation
191 {
192 template <class VectorType>
193 inline TmpVectors<VectorType>::TmpVectors(const unsigned int s_param,
195 : mem(vmem)
196 , data(s_param)
197 {}
198
199
200
201 template <class VectorType>
202 inline VectorType &TmpVectors<VectorType>::
203 operator[](const unsigned int i) const
204 {
205 AssertIndexRange(i, data.size());
206
207 Assert(data[i] != nullptr, ExcNotInitialized());
208 return *data[i];
209 }
210
211
212
213 template <class VectorType>
214 inline VectorType &
215 TmpVectors<VectorType>::operator()(const unsigned int i,
216 const VectorType & temp)
217 {
218 AssertIndexRange(i, data.size());
219 if (data[i] == nullptr)
220 {
221 data[i] = std::move(typename VectorMemory<VectorType>::Pointer(mem));
222 data[i]->reinit(temp);
223 }
224 return *data[i];
225 }
226 } // namespace SolverIDRImplementation
227} // namespace internal
228
229
230
231template <class VectorType>
234 const AdditionalData & data)
235 : SolverBase<VectorType>(cn, mem)
236 , additional_data(data)
237{}
238
239
240
241template <class VectorType>
242SolverIDR<VectorType>::SolverIDR(SolverControl &cn, const AdditionalData &data)
243 : SolverBase<VectorType>(cn)
244 , additional_data(data)
245{}
246
247
248
249template <class VectorType>
250void
251SolverIDR<VectorType>::print_vectors(const unsigned int,
252 const VectorType &,
253 const VectorType &,
254 const VectorType &) const
255{}
256
257
258
259template <class VectorType>
260template <typename MatrixType, typename PreconditionerType>
261void
262SolverIDR<VectorType>::solve(const MatrixType & A,
263 VectorType & x,
264 const VectorType & b,
265 const PreconditionerType &preconditioner)
266{
267 LogStream::Prefix prefix("IDR(s)");
268
270 unsigned int step = 0;
271
272 const unsigned int s = additional_data.s;
273
274 // Define temporary vectors which do not do not depend on s
275 typename VectorMemory<VectorType>::Pointer r_pointer(this->memory);
276 typename VectorMemory<VectorType>::Pointer v_pointer(this->memory);
277 typename VectorMemory<VectorType>::Pointer vhat_pointer(this->memory);
278 typename VectorMemory<VectorType>::Pointer uhat_pointer(this->memory);
279 typename VectorMemory<VectorType>::Pointer ghat_pointer(this->memory);
280
281 VectorType &r = *r_pointer;
282 VectorType &v = *v_pointer;
283 VectorType &vhat = *vhat_pointer;
284 VectorType &uhat = *uhat_pointer;
285 VectorType &ghat = *ghat_pointer;
286
287 r.reinit(x, true);
288 v.reinit(x, true);
289 vhat.reinit(x, true);
290 uhat.reinit(x, true);
291 ghat.reinit(x, true);
292
293 // Initial residual
294 A.vmult(r, x);
295 r.sadd(-1.0, 1.0, b);
296
297 // Check for convergent initial guess
298 double res = r.l2_norm();
299 iteration_state = this->iteration_status(step, res, x);
300 if (iteration_state == SolverControl::success)
301 return;
302
303 // Initialize sets of vectors/matrices whose size dependent on s
307 FullMatrix<double> M(s, s);
308
309 // Random number generator for vector entries of
310 // Q (normal distribution, mean=0 sigma=1)
311 std::mt19937 rng;
312 std::normal_distribution<> normal_distribution(0.0, 1.0);
313 for (unsigned int i = 0; i < s; ++i)
314 {
315 VectorType &tmp_g = G(i, x);
316 VectorType &tmp_u = U(i, x);
317 tmp_g = 0;
318 tmp_u = 0;
319
320 // Compute random set of s orthonormalized vectors Q
321 // Note: the first vector is chosen to be the initial
322 // residual to match BiCGStab (as is done in comparisons
323 // with BiCGStab in the papers listed in the documentation
324 // of this function)
325 VectorType &tmp_q = Q(i, x);
326 if (i != 0)
327 {
328 for (auto indx : tmp_q.locally_owned_elements())
329 tmp_q(indx) = normal_distribution(rng);
330 tmp_q.compress(VectorOperation::insert);
331 }
332 else
333 tmp_q = r;
334
335 for (unsigned int j = 0; j < i; ++j)
336 {
337 v = Q[j];
338 v *= (v * tmp_q) / (tmp_q * tmp_q);
339 tmp_q.add(-1.0, v);
340 }
341
342 if (i != 0)
343 tmp_q *= 1.0 / tmp_q.l2_norm();
344
345 M(i, i) = 1.;
346 }
347
348 double omega = 1.;
349
350 bool early_exit = false;
351
352 // Outer iteration
353 while (iteration_state == SolverControl::iterate)
354 {
355 ++step;
356
357 // Compute phi
358 Vector<double> phi(s);
359 for (unsigned int i = 0; i < s; ++i)
360 phi(i) = Q[i] * r;
361
362 // Inner iteration over s
363 for (unsigned int k = 0; k < s; ++k)
364 {
365 // Solve M(k:s)*gamma = phi(k:s)
366 Vector<double> gamma(s - k);
367 {
368 Vector<double> phik(s - k);
369 FullMatrix<double> Mk(s - k, s - k);
370 std::vector<unsigned int> indices;
371 unsigned int j = 0;
372 for (unsigned int i = k; i < s; ++i, ++j)
373 {
374 indices.push_back(i);
375 phik(j) = phi(i);
376 }
377 Mk.extract_submatrix_from(M, indices, indices);
378
379 FullMatrix<double> Mk_inv(s - k, s - k);
380 Mk_inv.invert(Mk);
381 Mk_inv.vmult(gamma, phik);
382 }
383
384 {
385 v = r;
386
387 unsigned int j = 0;
388 for (unsigned int i = k; i < s; ++i, ++j)
389 v.add(-1.0 * gamma(j), G[i]);
390 preconditioner.vmult(vhat, v);
391
392 uhat = vhat;
393 uhat *= omega;
394 j = 0;
395 for (unsigned int i = k; i < s; ++i, ++j)
396 uhat.add(gamma(j), U[i]);
397 A.vmult(ghat, uhat);
398 }
399
400 // Update G and U
401 // Orthogonalize ghat to Q0,..,Q_{k-1}
402 // and update uhat
403 for (unsigned int i = 0; i < k; ++i)
404 {
405 double alpha = (Q[i] * ghat) / M(i, i);
406 ghat.add(-alpha, G[i]);
407 uhat.add(-alpha, U[i]);
408 }
409 G[k] = ghat;
410 U[k] = uhat;
411
412 // Update kth column of M
413 for (unsigned int i = k; i < s; ++i)
414 M(i, k) = Q[i] * G[k];
415
416 // Orthogonalize r to Q0,...,Qk,
417 // update x
418 {
419 double beta = phi(k) / M(k, k);
420 r.add(-1.0 * beta, G[k]);
421 x.add(beta, U[k]);
422
423 print_vectors(step, x, r, U[k]);
424
425 // Check for early convergence. If so, store
426 // information in early_exit so that outer iteration
427 // is broken before recomputing the residual
428 res = r.l2_norm();
429 iteration_state = this->iteration_status(step, res, x);
430 if (iteration_state != SolverControl::iterate)
431 {
432 early_exit = true;
433 break;
434 }
435
436 // Update phi
437 if (k + 1 < s)
438 {
439 for (unsigned int i = 0; i < k + 1; ++i)
440 phi(i) = 0.0;
441 for (unsigned int i = k + 1; i < s; ++i)
442 phi(i) -= beta * M(i, k);
443 }
444 }
445 }
446 if (early_exit == true)
447 break;
448
449 // Update r and x
450 preconditioner.vmult(vhat, r);
451 A.vmult(v, vhat);
452
453 omega = (v * r) / (v * v);
454
455 r.add(-1.0 * omega, v);
456 x.add(omega, vhat);
457
458 print_vectors(step, x, r, vhat);
459
460 // Check for convergence
461 res = r.l2_norm();
462 iteration_state = this->iteration_status(step, res, x);
463 if (iteration_state != SolverControl::iterate)
464 break;
465 }
466
467 if (iteration_state != SolverControl::success)
468 AssertThrow(false, SolverControl::NoConvergence(step, res));
469}
470
471
472#endif // DOXYGEN
473
475
476#endif
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
SolverIDR(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverIDR(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
AdditionalData additional_data
Definition: solver_idr.h:179
virtual ~SolverIDR() override=default
Definition: vector.h:110
TmpVectors(const unsigned int s_param, VectorMemory< VectorType > &vmem)
VectorType & operator()(const unsigned int i, const VectorType &temp)
VectorType & operator[](const unsigned int i) const
VectorMemory< VectorType > & mem
Definition: solver_idr.h:83
std::vector< typename VectorMemory< VectorType >::Pointer > data
Definition: solver_idr.h:88
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
static const char U
static const char A
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)
const unsigned int s
Definition: solver_idr.h:131
AdditionalData(const unsigned int s=2)
Definition: solver_idr.h:127