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