Dune Core Modules (2.9.1)

multitypeblockmatrix.hh
1// SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_ISTL_MULTITYPEBLOCKMATRIX_HH
6#define DUNE_ISTL_MULTITYPEBLOCKMATRIX_HH
7
8#include <cmath>
9#include <iostream>
10#include <tuple>
11
12#include <dune/common/hybridutilities.hh>
13
14#include "istlexception.hh"
15
16// forward declaration
17namespace Dune
18{
19 template<typename FirstRow, typename... Args>
20 class MultiTypeBlockMatrix;
21
22 template<int I, int crow, int remain_row>
23 class MultiTypeBlockMatrix_Solver;
24}
25
26#include "gsetc.hh"
27
28namespace Dune {
29
43 template<typename FirstRow, typename... Args>
45 : public std::tuple<FirstRow, Args...>
46 {
47 using ParentType = std::tuple<FirstRow, Args...>;
48 public:
49
51 using ParentType::ParentType;
52
56 typedef MultiTypeBlockMatrix<FirstRow, Args...> type;
57
59 using size_type = std::size_t;
60
61 typedef typename FirstRow::field_type field_type;
62
64 static constexpr size_type N()
65 {
66 return 1+sizeof...(Args);
67 }
68
74 [[deprecated("Use method 'N' instead")]]
75 static constexpr size_type size()
76 {
77 return 1+sizeof...(Args);
78 }
79
81 static constexpr size_type M()
82 {
83 return FirstRow::size();
84 }
85
102 template< size_type index >
103 auto
104 operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable)
105 -> decltype(std::get<index>(*this))
106 {
107 return std::get<index>(*this);
108 }
109
115 template< size_type index >
116 auto
117 operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable) const
118 -> decltype(std::get<index>(*this))
119 {
120 return std::get<index>(*this);
121 }
122
126 template<typename T>
127 void operator= (const T& newval) {
128 using namespace Dune::Hybrid;
129 auto size = index_constant<1+sizeof...(Args)>();
130 // Since Dune::Hybrid::size(MultiTypeBlockMatrix) is not implemented,
131 // we cannot use a plain forEach(*this, ...). This could be achieved,
132 // e.g., by implementing a static size() method.
133 forEach(integralRange(size), [&](auto&& i) {
134 (*this)[i] = newval;
135 });
136 }
137
138 //===== vector space arithmetic
139
141 MultiTypeBlockMatrix& operator*= (const field_type& k)
142 {
143 auto size = index_constant<N()>();
145 (*this)[i] *= k;
146 });
147
148 return *this;
149 }
150
152 MultiTypeBlockMatrix& operator/= (const field_type& k)
153 {
154 auto size = index_constant<N()>();
156 (*this)[i] /= k;
157 });
158
159 return *this;
160 }
161
162
169 {
170 auto size = index_constant<N()>();
172 (*this)[i] += b[i];
173 });
174
175 return *this;
176 }
177
184 {
185 auto size = index_constant<N()>();
187 (*this)[i] -= b[i];
188 });
189
190 return *this;
191 }
192
195 template<typename X, typename Y>
196 void mv (const X& x, Y& y) const {
197 static_assert(X::size() == M(), "length of x does not match row length");
198 static_assert(Y::size() == N(), "length of y does not match row count");
199 y = 0; //reset y (for mv uses umv)
200 umv(x,y);
201 }
202
205 template<typename X, typename Y>
206 void umv (const X& x, Y& y) const {
207 static_assert(X::size() == M(), "length of x does not match row length");
208 static_assert(Y::size() == N(), "length of y does not match row count");
209 using namespace Dune::Hybrid;
210 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
211 using namespace Dune::Hybrid; // needed for icc, see issue #31
212 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
213 (*this)[i][j].umv(x[j], y[i]);
214 });
215 });
216 }
217
220 template<typename X, typename Y>
221 void mmv (const X& x, Y& y) const {
222 static_assert(X::size() == M(), "length of x does not match row length");
223 static_assert(Y::size() == N(), "length of y does not match row count");
224 using namespace Dune::Hybrid;
225 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
226 using namespace Dune::Hybrid; // needed for icc, see issue #31
227 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
228 (*this)[i][j].mmv(x[j], y[i]);
229 });
230 });
231 }
232
235 template<typename AlphaType, typename X, typename Y>
236 void usmv (const AlphaType& alpha, const X& x, Y& y) const {
237 static_assert(X::size() == M(), "length of x does not match row length");
238 static_assert(Y::size() == N(), "length of y does not match row count");
239 using namespace Dune::Hybrid;
240 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
241 using namespace Dune::Hybrid; // needed for icc, see issue #31
242 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
243 (*this)[i][j].usmv(alpha, x[j], y[i]);
244 });
245 });
246 }
247
250 template<typename X, typename Y>
251 void mtv (const X& x, Y& y) const {
252 static_assert(X::size() == N(), "length of x does not match number of rows");
253 static_assert(Y::size() == M(), "length of y does not match number of columns");
254 y = 0;
255 umtv(x,y);
256 }
257
260 template<typename X, typename Y>
261 void umtv (const X& x, Y& y) const {
262 static_assert(X::size() == N(), "length of x does not match number of rows");
263 static_assert(Y::size() == M(), "length of y does not match number of columns");
264 using namespace Dune::Hybrid;
265 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
266 using namespace Dune::Hybrid; // needed for icc, see issue #31
267 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
268 (*this)[j][i].umtv(x[j], y[i]);
269 });
270 });
271 }
272
275 template<typename X, typename Y>
276 void mmtv (const X& x, Y& y) const {
277 static_assert(X::size() == N(), "length of x does not match number of rows");
278 static_assert(Y::size() == M(), "length of y does not match number of columns");
279 using namespace Dune::Hybrid;
280 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
281 using namespace Dune::Hybrid; // needed for icc, see issue #31
282 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
283 (*this)[j][i].mmtv(x[j], y[i]);
284 });
285 });
286 }
287
290 template<typename X, typename Y>
291 void usmtv (const field_type& alpha, const X& x, Y& y) const {
292 static_assert(X::size() == N(), "length of x does not match number of rows");
293 static_assert(Y::size() == M(), "length of y does not match number of columns");
294 using namespace Dune::Hybrid;
295 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
296 using namespace Dune::Hybrid; // needed for icc, see issue #31
297 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
298 (*this)[j][i].usmtv(alpha, x[j], y[i]);
299 });
300 });
301 }
302
305 template<typename X, typename Y>
306 void umhv (const X& x, Y& y) const {
307 static_assert(X::size() == N(), "length of x does not match number of rows");
308 static_assert(Y::size() == M(), "length of y does not match number of columns");
309 using namespace Dune::Hybrid;
310 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
311 using namespace Dune::Hybrid; // needed for icc, see issue #31
312 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
313 (*this)[j][i].umhv(x[j], y[i]);
314 });
315 });
316 }
317
320 template<typename X, typename Y>
321 void mmhv (const X& x, Y& y) const {
322 static_assert(X::size() == N(), "length of x does not match number of rows");
323 static_assert(Y::size() == M(), "length of y does not match number of columns");
324 using namespace Dune::Hybrid;
325 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
326 using namespace Dune::Hybrid; // needed for icc, see issue #31
327 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
328 (*this)[j][i].mmhv(x[j], y[i]);
329 });
330 });
331 }
332
335 template<typename X, typename Y>
336 void usmhv (const field_type& alpha, const X& x, Y& y) const {
337 static_assert(X::size() == N(), "length of x does not match number of rows");
338 static_assert(Y::size() == M(), "length of y does not match number of columns");
339 using namespace Dune::Hybrid;
340 forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
341 using namespace Dune::Hybrid; // needed for icc, see issue #31
342 forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
343 (*this)[j][i].usmhv(alpha, x[j], y[i]);
344 });
345 });
346 }
347
348
349 //===== norms
350
352 auto frobenius_norm2 () const
353 {
354 using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
355 typename FieldTraits<field_type_00>::real_type sum=0;
356
357 auto rows = index_constant<N()>();
358 Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
360 Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
361 sum += (*this)[i][j].frobenius_norm2();
362 });
363 });
364
365 return sum;
366 }
367
369 typename FieldTraits<field_type>::real_type frobenius_norm () const
370 {
371 return sqrt(frobenius_norm2());
372 }
373
375 auto infinity_norm () const
376 {
377 using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
378 using std::max;
379 typename FieldTraits<field_type_00>::real_type norm=0;
380
381 auto rows = index_constant<N()>();
382 Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
383
384 typename FieldTraits<field_type_00>::real_type sum=0;
386 Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
387 sum += (*this)[i][j].infinity_norm();
388 });
389 norm = max(sum, norm);
390 });
391
392 return norm;
393 }
394
396 auto infinity_norm_real () const
397 {
398 using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
399 using std::max;
400 typename FieldTraits<field_type_00>::real_type norm=0;
401
402 auto rows = index_constant<N()>();
403 Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
404
405 typename FieldTraits<field_type_00>::real_type sum=0;
407 Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
408 sum += (*this)[i][j].infinity_norm_real();
409 });
410 norm = max(sum, norm);
411 });
412
413 return norm;
414 }
415
416 };
417
423 template<typename T1, typename... Args>
424 std::ostream& operator<< (std::ostream& s, const MultiTypeBlockMatrix<T1,Args...>& m) {
427 using namespace Dune::Hybrid;
428 forEach(integralRange(N), [&](auto&& i) {
429 using namespace Dune::Hybrid; // needed for icc, see issue #31
430 forEach(integralRange(M), [&](auto&& j) {
431 s << "\t(" << i << ", " << j << "): \n" << m[i][j];
432 });
433 });
434 s << std::endl;
435 return s;
436 }
437
438 //make algmeta_itsteps known
439 template<int I, typename M>
440 struct algmeta_itsteps;
441
448 template<int I, int crow, int ccol, int remain_col> //MultiTypeBlockMatrix_Solver_Col: iterating over one row
449 class MultiTypeBlockMatrix_Solver_Col { //calculating b- A[i][j]*x[j]
450 public:
454 template <typename Trhs, typename TVector, typename TMatrix, typename K>
455 static void calc_rhs(const TMatrix& A, TVector& x, TVector& v, Trhs& b, const K& w) {
456 std::get<ccol>( std::get<crow>(A) ).mmv( std::get<ccol>(x), b );
458 }
459
460 };
461 template<int I, int crow, int ccol> //MultiTypeBlockMatrix_Solver_Col recursion end
462 class MultiTypeBlockMatrix_Solver_Col<I,crow,ccol,0> {
463 public:
464 template <typename Trhs, typename TVector, typename TMatrix, typename K>
465 static void calc_rhs(const TMatrix&, TVector&, TVector&, Trhs&, const K&) {}
466 };
467
468
469
476 template<int I, int crow, int remain_row>
478 public:
479
483 template <typename TVector, typename TMatrix, typename K>
484 static void dbgs(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
485 TVector xold(x);
486 xold=x; //store old x values
488 x *= w;
489 x.axpy(1-w,xold); //improve x
490 }
491 template <typename TVector, typename TMatrix, typename K>
492 static void dbgs(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
493 auto rhs = std::get<crow> (b);
494
495 MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
496 //solve on blocklevel I-1
497 using M =
498 typename std::remove_cv<
499 typename std::remove_reference<
500 decltype(std::get<crow>( std::get<crow>(A)))
501 >::type
502 >::type;
503 algmeta_itsteps<I-1,M>::dbgs(std::get<crow>( std::get<crow>(A)), std::get<crow>(x),rhs,w);
505 }
506
507
508
512 template <typename TVector, typename TMatrix, typename K>
513 static void bsorf(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
514 TVector v;
515 v=x; //use latest x values in right side calculation
517
518 }
519 template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
520 static void bsorf(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
521 auto rhs = std::get<crow> (b);
522
523 MultiTypeBlockMatrix_Solver_Col<I,crow,0,TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
524 //solve on blocklevel I-1
525 using M =
526 typename std::remove_cv<
527 typename std::remove_reference<
528 decltype(std::get<crow>( std::get<crow>(A)))
529 >::type
530 >::type;
531 algmeta_itsteps<I-1,M>::bsorf(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
532 std::get<crow>(x).axpy(w,std::get<crow>(v));
534 }
535
539 template <typename TVector, typename TMatrix, typename K>
540 static void bsorb(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
541 TVector v;
542 v=x; //use latest x values in right side calculation
544
545 }
546 template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
547 static void bsorb(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
548 auto rhs = std::get<crow> (b);
549
550 MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
551 //solve on blocklevel I-1
552 using M =
553 typename std::remove_cv<
554 typename std::remove_reference<
555 decltype(std::get<crow>( std::get<crow>(A)))
556 >::type
557 >::type;
558 algmeta_itsteps<I-1,M>::bsorb(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
559 std::get<crow>(x).axpy(w,std::get<crow>(v));
561 }
562
563
567 template <typename TVector, typename TMatrix, typename K>
568 static void dbjac(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
569 TVector v(x);
570 v=0; //calc new x in v
572 x.axpy(w,v); //improve x
573 }
574 template <typename TVector, typename TMatrix, typename K>
575 static void dbjac(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
576 auto rhs = std::get<crow> (b);
577
578 MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
579 //solve on blocklevel I-1
580 using M =
581 typename std::remove_cv<
582 typename std::remove_reference<
583 decltype(std::get<crow>( std::get<crow>(A)))
584 >::type
585 >::type;
586 algmeta_itsteps<I-1,M>::dbjac(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
588 }
589
590
591
592
593 };
594 template<int I, int crow> //recursion end for remain_row = 0
595 class MultiTypeBlockMatrix_Solver<I,crow,0> {
596 public:
597 template <typename TVector, typename TMatrix, typename K>
598 static void dbgs(const TMatrix&, TVector&, TVector&,
599 const TVector&, const K&) {}
600
601 template <typename TVector, typename TMatrix, typename K>
602 static void bsorf(const TMatrix&, TVector&, TVector&,
603 const TVector&, const K&) {}
604
605 template <typename TVector, typename TMatrix, typename K>
606 static void bsorb(const TMatrix&, TVector&, TVector&,
607 const TVector&, const K&) {}
608
609 template <typename TVector, typename TMatrix, typename K>
610 static void dbjac(const TMatrix&, TVector&, TVector&,
611 const TVector&, const K&) {}
612 };
613
614} // end namespace Dune
615
616namespace std
617{
622 template <size_t i, typename... Args>
623 struct tuple_element<i,Dune::MultiTypeBlockMatrix<Args...> >
624 {
625 using type = typename std::tuple_element<i, std::tuple<Args...> >::type;
626 };
627}
628#endif
part of solvers for MultiTypeBlockVector & MultiTypeBlockMatrix types
Definition: multitypeblockmatrix.hh:449
solver for MultiTypeBlockVector & MultiTypeBlockMatrix types
Definition: multitypeblockmatrix.hh:477
A Matrix class to support different block types.
Definition: multitypeblockmatrix.hh:46
constexpr index_constant< 0 > _0
Compile time index with value 0.
Definition: indices.hh:53
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:30
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:268
constexpr auto integralRange(const Begin &begin, const End &end)
Create an integral range.
Definition: hybridutilities.hh:184
void bsorb(const M &A, X &x, const Y &b, const K &w)
SSOR step.
Definition: gsetc.hh:646
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:658
void dbgs(const M &A, X &x, const Y &b, const K &w)
GS step.
Definition: gsetc.hh:622
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:634
MultiTypeBlockMatrix< FirstRow, Args... > type
Definition: multitypeblockmatrix.hh:56
static void dbjac(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:568
MultiTypeBlockMatrix & operator+=(const MultiTypeBlockMatrix &b)
Add the entries of another matrix to this one.
Definition: multitypeblockmatrix.hh:168
void mv(const X &x, Y &y) const
y = A x
Definition: multitypeblockmatrix.hh:196
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: multitypeblockmatrix.hh:276
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: multitypeblockmatrix.hh:321
static void dbgs(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:484
static void bsorb(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:540
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: multitypeblockmatrix.hh:336
void usmv(const AlphaType &alpha, const X &x, Y &y) const
y += alpha A x
Definition: multitypeblockmatrix.hh:236
MultiTypeBlockMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: multitypeblockmatrix.hh:152
void umhv(const X &x, Y &y) const
y += A^H x
Definition: multitypeblockmatrix.hh:306
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: multitypeblockmatrix.hh:369
static constexpr size_type N()
Return the number of matrix rows.
Definition: multitypeblockmatrix.hh:64
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: multitypeblockmatrix.hh:291
void operator=(const T &newval)
Definition: multitypeblockmatrix.hh:127
void umv(const X &x, Y &y) const
y += A x
Definition: multitypeblockmatrix.hh:206
static void calc_rhs(const TMatrix &A, TVector &x, TVector &v, Trhs &b, const K &w)
Definition: multitypeblockmatrix.hh:455
static void bsorf(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:513
static constexpr size_type size()
Return the number of matrix rows.
Definition: multitypeblockmatrix.hh:75
auto infinity_norm_real() const
Bastardized version of the infinity-norm / row-sum norm.
Definition: multitypeblockmatrix.hh:396
auto frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: multitypeblockmatrix.hh:352
std::size_t size_type
Type used for sizes.
Definition: multitypeblockmatrix.hh:59
auto infinity_norm() const
Bastardized version of the infinity-norm / row-sum norm.
Definition: multitypeblockmatrix.hh:375
MultiTypeBlockMatrix & operator-=(const MultiTypeBlockMatrix &b)
Subtract the entries of another matrix from this one.
Definition: multitypeblockmatrix.hh:183
auto operator[](const std::integral_constant< size_type, index > indexVariable) -> decltype(std::get< index >(*this))
Random-access operator.
Definition: multitypeblockmatrix.hh:104
void mtv(const X &x, Y &y) const
y = A^T x
Definition: multitypeblockmatrix.hh:251
static constexpr size_type M()
Return the number of matrix columns.
Definition: multitypeblockmatrix.hh:81
void mmv(const X &x, Y &y) const
y -= A x
Definition: multitypeblockmatrix.hh:221
void umtv(const X &x, Y &y) const
y += A^T x
Definition: multitypeblockmatrix.hh:261
MultiTypeBlockMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: multitypeblockmatrix.hh:141
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Dune namespace.
Definition: alignedallocator.hh:13
STL namespace.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)