version 3.10-dev
parallelmatrixadapter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_LINEAR_PARALLEL_MATRIX_ADAPTER_HH
13#define DUMUX_LINEAR_PARALLEL_MATRIX_ADAPTER_HH
14
15#include <dune/common/hybridutilities.hh>
16#include <dune/istl/operators.hh>
17#include <dune/istl/bcrsmatrix.hh>
19
20namespace Dumux {
21
26template<class M, class X, class Y>
27class ParallelMultiTypeMatrixAdapter : public Dune::MatrixAdapter<M,X,Y>
28{
29 using ParentType = Dune::MatrixAdapter<M,X,Y>;
30public:
32 typedef M matrix_type;
33 typedef X domain_type;
34 typedef Y range_type;
35 typedef typename X::field_type field_type;
36
38 explicit ParallelMultiTypeMatrixAdapter (const M& A)
39 : ParentType(A) {}
40
42 explicit ParallelMultiTypeMatrixAdapter (std::shared_ptr<const M> A)
43 : ParentType(A) {}
44
46 void apply (const X& x, Y& y) const override
47 {
48 y = 0.0;
49
50 auto& A = this->getmat();
51 using namespace Dune::Hybrid;
52 forEach(integralRange(Dune::Hybrid::size(y)), [&](auto&& i)
53 {
54 forEach(integralRange(Dune::Hybrid::size(x)), [&](auto&& j)
55 {
56 const auto& mat = A[i][j];
57 const auto& xx = x[j];
58 auto& yy = y[i];
59
60 Dumux::parallelFor(mat.N(), [&](const std::size_t ii)
61 {
62 const auto& row = mat[ii];
63 const auto endj = row.end();
64 for (auto jj=row.begin(); jj!=endj; ++jj)
65 {
66 const auto& xj = Dune::Impl::asVector(xx[jj.index()]);
67 auto&& yi = Dune::Impl::asVector(yy[ii]);
68 Dune::Impl::asMatrix(*jj).umv(xj, yi);
69 }
70 });
71 });
72 });
73 }
74
76 void applyscaleadd (field_type alpha, const X& x, Y& y) const override
77 {
78 auto& A = this->getmat();
79 using namespace Dune::Hybrid;
80 forEach(integralRange(Dune::Hybrid::size(y)), [&](auto&& i)
81 {
82 forEach(integralRange(Dune::Hybrid::size(x)), [&](auto&& j)
83 {
84 const auto& mat = A[i][j];
85 const auto& xx = x[j];
86 auto& yy = y[i];
87
88 Dumux::parallelFor(mat.N(), [&](const std::size_t ii)
89 {
90 const auto& row = mat[ii];
91 const auto endj = row.end();
92 for (auto jj=row.begin(); jj!=endj; ++jj)
93 {
94 const auto& xj = Dune::Impl::asVector(xx[jj.index()]);
95 auto&& yi = Dune::Impl::asVector(yy[ii]);
96 Dune::Impl::asMatrix(*jj).usmv(alpha, xj, yi);
97 }
98 });
99 });
100 });
101 }
102};
103
104} // end namespace Dumux
105
106#endif
Adapter to turn a multi-type matrix into a thread-parallel linear operator. Adapts a matrix to the as...
Definition: parallelmatrixadapter.hh:28
X domain_type
Definition: parallelmatrixadapter.hh:33
X::field_type field_type
Definition: parallelmatrixadapter.hh:35
ParallelMultiTypeMatrixAdapter(const M &A)
constructor: just store a reference to a matrix
Definition: parallelmatrixadapter.hh:38
void apply(const X &x, Y &y) const override
apply operator to x:
Definition: parallelmatrixadapter.hh:46
Y range_type
Definition: parallelmatrixadapter.hh:34
ParallelMultiTypeMatrixAdapter(std::shared_ptr< const M > A)
constructor: store an std::shared_ptr to a matrix
Definition: parallelmatrixadapter.hh:42
M matrix_type
export types
Definition: parallelmatrixadapter.hh:32
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
apply operator to x, scale and add:
Definition: parallelmatrixadapter.hh:76
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
Definition: adapt.hh:17
Parallel for loop (multithreading)