deal.II version 9.7.0
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_solver_idr_h
16#define dealii_solver_idr_h
17
18
19#include <deal.II/base/config.h>
20
25
28#include <deal.II/lac/solver.h>
30
31#include <cmath>
32#include <random>
33
35
40
41namespace internal
42{
47 {
52 template <typename VectorType>
54 {
55 public:
59 TmpVectors(const unsigned int s_param, VectorMemory<VectorType> &vmem);
60
64 ~TmpVectors() = default;
65
70 VectorType &
71 operator[](const unsigned int i) const;
72
79 VectorType &
80 operator()(const unsigned int i, const VectorType &temp);
81
82 private:
87
91 std::vector<typename VectorMemory<VectorType>::Pointer> data;
92 };
93 } // namespace SolverIDRImplementation
94} // namespace internal
95
118template <typename VectorType = Vector<double>>
120class SolverIDR : public SolverBase<VectorType>
121{
122public:
127 {
131 explicit AdditionalData(const unsigned int s = 2)
132 : s(s)
133 {}
134
135 const unsigned int s;
136 };
137
143 const AdditionalData &data = AdditionalData());
144
150 const AdditionalData &data = AdditionalData());
151
155 virtual ~SolverIDR() override = default;
156
160 template <typename MatrixType, typename PreconditionerType>
164 void solve(const MatrixType &A,
165 VectorType &x,
166 const VectorType &b,
167 const PreconditionerType &preconditioner);
168
169protected:
175 virtual void
176 print_vectors(const unsigned int step,
177 const VectorType &x,
178 const VectorType &r,
179 const VectorType &d) const;
180
181private:
186};
187
189/*------------------------- Implementation ----------------------------*/
190
191#ifndef DOXYGEN
192
193
194namespace internal
195{
196 namespace SolverIDRImplementation
197 {
198 template <typename VectorType>
199 inline TmpVectors<VectorType>::TmpVectors(const unsigned int s_param,
201 : mem(vmem)
202 , data(s_param)
203 {}
204
205
206
207 template <typename VectorType>
208 inline VectorType &
209 TmpVectors<VectorType>::operator[](const unsigned int i) const
210 {
211 AssertIndexRange(i, data.size());
212
213 Assert(data[i] != nullptr, ExcNotInitialized());
214 return *data[i];
215 }
216
217
218
219 template <typename VectorType>
220 inline VectorType &
221 TmpVectors<VectorType>::operator()(const unsigned int i,
222 const VectorType &temp)
223 {
224 AssertIndexRange(i, data.size());
225 if (data[i] == nullptr)
226 {
227 data[i] = std::move(typename VectorMemory<VectorType>::Pointer(mem));
228 data[i]->reinit(temp, true);
229 }
230 return *data[i];
231 }
232
233
234
235 template <typename VectorType,
236 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
237 * = nullptr>
238 unsigned int
239 n_blocks(const VectorType &)
240 {
241 return 1;
242 }
243
244
245
246 template <typename VectorType,
247 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
248 nullptr>
249 unsigned int
250 n_blocks(const VectorType &vector)
251 {
252 return vector.n_blocks();
253 }
254
255
256
257 template <typename VectorType,
258 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
259 * = nullptr>
260 VectorType &
261 block(VectorType &vector, const unsigned int b)
262 {
263 AssertDimension(b, 0);
264 return vector;
265 }
266
267
268
269 template <typename VectorType,
270 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
271 nullptr>
272 typename VectorType::BlockType &
273 block(VectorType &vector, const unsigned int b)
274 {
275 return vector.block(b);
276 }
277
278 } // namespace SolverIDRImplementation
279} // namespace internal
280
281
282
283template <typename VectorType>
287 const AdditionalData &data)
288 : SolverBase<VectorType>(cn, mem)
289 , additional_data(data)
290{}
291
292
293
294template <typename VectorType>
296SolverIDR<VectorType>::SolverIDR(SolverControl &cn, const AdditionalData &data)
297 : SolverBase<VectorType>(cn)
298 , additional_data(data)
299{}
300
301
302
303template <typename VectorType>
305void SolverIDR<VectorType>::print_vectors(const unsigned int,
306 const VectorType &,
307 const VectorType &,
308 const VectorType &) const
309{}
310
311
312
313template <typename VectorType>
315template <typename MatrixType, typename PreconditionerType>
319void SolverIDR<VectorType>::solve(const MatrixType &A,
320 VectorType &x,
321 const VectorType &b,
322 const PreconditionerType &preconditioner)
323{
324 LogStream::Prefix prefix("IDR(s)");
325
327 unsigned int step = 0;
328
329 const unsigned int s = additional_data.s;
330
331 // Define temporary vectors which do not do not depend on s
332 typename VectorMemory<VectorType>::Pointer r_pointer(this->memory);
333 typename VectorMemory<VectorType>::Pointer v_pointer(this->memory);
334 typename VectorMemory<VectorType>::Pointer uhat_pointer(this->memory);
335
336 VectorType &r = *r_pointer;
337 VectorType &v = *v_pointer;
338 VectorType &uhat = *uhat_pointer;
339
340 r.reinit(x, true);
341 v.reinit(x, true);
342 uhat.reinit(x, true);
343
344 // Initial residual
345 A.vmult(r, x);
346 r.sadd(-1.0, 1.0, b);
347
348 using value_type = typename VectorType::value_type;
349 using real_type = typename numbers::NumberTraits<value_type>::real_type;
350
351 // Check for convergent initial guess
352 real_type res = r.l2_norm();
353 iteration_state = this->iteration_status(step, res, x);
354 if (iteration_state == SolverControl::success)
355 return;
356
357 // Initialize sets of vectors/matrices whose size dependent on s
362
363 // Random number generator for vector entries of
364 // Q (normal distribution, mean=0 sigma=1)
365 std::mt19937 rng;
366 std::normal_distribution<> normal_distribution(0.0, 1.0);
367 for (unsigned int i = 0; i < s; ++i)
368 {
369 // Initialize vectors
370 G(i, x);
371 U(i, x);
372
373 // Compute random set of s orthonormalized vectors Q
374 // Note: the first vector is chosen to be the initial
375 // residual to match BiCGStab (as is done in comparisons
376 // with BiCGStab in the papers listed in the documentation
377 // of this function)
378 VectorType &tmp_q = Q(i, x);
379 if (i != 0)
380 {
381 for (unsigned int b = 0;
382 b < internal::SolverIDRImplementation::n_blocks(tmp_q);
383 ++b)
384 for (auto indx : internal::SolverIDRImplementation::block(tmp_q, b)
385 .locally_owned_elements())
386 internal::SolverIDRImplementation::block(tmp_q, b)(indx) =
387 normal_distribution(rng);
388 tmp_q.compress(VectorOperation::insert);
389 }
390 else
391 tmp_q = r;
392
393 for (unsigned int j = 0; j < i; ++j)
394 {
395 v = Q[j];
396 v *= (v * tmp_q) / (tmp_q * tmp_q);
397 tmp_q.add(-1.0, v);
398 }
399
400 if (i != 0)
401 tmp_q *= 1.0 / tmp_q.l2_norm();
402
403 M(i, i) = 1.;
404 }
405
406 value_type omega = 1.;
407
408 bool early_exit = false;
409
410 // Outer iteration
411 while (iteration_state == SolverControl::iterate)
412 {
413 ++step;
414
415 // Compute phi
416 Vector<value_type> phi(s);
417 for (unsigned int i = 0; i < s; ++i)
418 phi(i) = Q[i] * r;
419
420 // Inner iteration over s
421 for (unsigned int k = 0; k < s; ++k)
422 {
423 // Solve M(k:s)*gamma = phi(k:s)
425 {
426 Vector<value_type> phik(s - k);
427 FullMatrix<value_type> Mk(s - k, s - k);
428 std::vector<unsigned int> indices;
429 unsigned int j = 0;
430 for (unsigned int i = k; i < s; ++i, ++j)
431 {
432 indices.push_back(i);
433 phik(j) = phi(i);
434 }
435 Mk.extract_submatrix_from(M, indices, indices);
436
437 FullMatrix<value_type> Mk_inv(s - k, s - k);
438 Mk_inv.invert(Mk);
439 Mk_inv.vmult(gamma, phik);
440 }
441
442 v = r;
443
444 if (step > 1)
445 {
446 for (unsigned int i = k, j = 0; i < s; ++i, ++j)
447 v.add(-gamma(j), G[i]);
448 }
449
450 preconditioner.vmult(uhat, v);
451
452 if (step > 1)
453 {
454 uhat.sadd(omega, gamma(0), U[k]);
455 for (unsigned int i = k + 1, j = 1; i < s; ++i, ++j)
456 uhat.add(gamma(j), U[i]);
457 }
458 else
459 uhat *= omega;
460
461 A.vmult(G[k], uhat);
462
463 // Update G and U
464 // Orthogonalize G[k] to Q0,..,Q_{k-1} and update uhat
465 if (k > 0)
466 {
467 value_type alpha = Q[0] * G[k] / M(0, 0);
468 for (unsigned int i = 1; i < k; ++i)
469 {
470 const value_type alpha_old = alpha;
471 alpha = G[k].add_and_dot(-alpha, G[i - 1], Q[i]) / M(i, i);
472
473 // update uhat every other iteration to reduce vector access
474 if (i % 2 == 1)
475 uhat.add(-alpha_old, U[i - 1], -alpha, U[i]);
476 }
477 M(k, k) = G[k].add_and_dot(-alpha, G[k - 1], Q[k]);
478 if (k % 2 == 1)
479 uhat.add(-alpha, U[k - 1]);
480 }
481 else
482 M(k, k) = G[k] * Q[k];
483
484 U[k].swap(uhat);
485
486 // Update kth column of M
487 for (unsigned int i = k + 1; i < s; ++i)
488 M(i, k) = Q[i] * G[k];
489
490 // Orthogonalize r to Q0,...,Qk, update x
491 {
492 const value_type beta = phi(k) / M(k, k);
493 r.add(-beta, G[k]);
494 x.add(beta, U[k]);
495
496 print_vectors(step, x, r, U[k]);
497
498 // Check for early convergence. If so, store
499 // information in early_exit so that outer iteration
500 // is broken before recomputing the residual
501 res = r.l2_norm();
502 iteration_state = this->iteration_status(step, res, x);
503 if (iteration_state != SolverControl::iterate)
504 {
505 early_exit = true;
506 break;
507 }
508
509 // Update phi
510 if (k + 1 < s)
511 {
512 for (unsigned int i = 0; i < k + 1; ++i)
513 phi(i) = 0.0;
514 for (unsigned int i = k + 1; i < s; ++i)
515 phi(i) -= beta * M(i, k);
516 }
517 }
518 }
519 if (early_exit == true)
520 break;
521
522 // Update r and x
523 preconditioner.vmult(uhat, r);
524 A.vmult(v, uhat);
525
526 omega = (v * r) / (v * v);
527
528 res = std::sqrt(r.add_and_dot(-1.0 * omega, v, r));
529 x.add(omega, uhat);
530
531 print_vectors(step, x, r, uhat);
532
533 // Check for convergence
534 iteration_state = this->iteration_status(step, res, x);
535 if (iteration_state != SolverControl::iterate)
536 break;
537 }
538
539 if (iteration_state != SolverControl::success)
540 AssertThrow(false, SolverControl::NoConvergence(step, res));
541}
542
543
544#endif // DOXYGEN
545
547
548#endif
SolverBase(SolverControl &solver_control, VectorMemory< VectorType > &vector_memory)
@ 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:185
virtual ~SolverIDR() override=default
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
std::vector< typename VectorMemory< VectorType >::Pointer > data
Definition solver_idr.h:91
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:248
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
constexpr char U
constexpr char A
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
long double gamma(const unsigned int n)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
const unsigned int s
Definition solver_idr.h:135
AdditionalData(const unsigned int s=2)
Definition solver_idr.h:131