Line | Branch | Exec | Source |
---|---|---|---|
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Linear | ||
10 | * \brief A parallel version of a linear operator | ||
11 | */ | ||
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> | ||
18 | #include <dumux/parallel/parallel_for.hh> | ||
19 | |||
20 | namespace Dumux { | ||
21 | |||
22 | /*! | ||
23 | * \brief Adapter to turn a multi-type matrix into a thread-parallel linear operator. | ||
24 | * Adapts a matrix to the assembled linear operator interface | ||
25 | */ | ||
26 | template<class M, class X, class Y> | ||
27 |
1/4✓ Branch 0 taken 2844 times.
✗ Branch 1 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
5688 | class ParallelMultiTypeMatrixAdapter : public Dune::MatrixAdapter<M,X,Y> |
28 | { | ||
29 | using ParentType = Dune::MatrixAdapter<M,X,Y>; | ||
30 | public: | ||
31 | //! export types | ||
32 | typedef M matrix_type; | ||
33 | typedef X domain_type; | ||
34 | typedef Y range_type; | ||
35 | typedef typename X::field_type field_type; | ||
36 | |||
37 | //! constructor: just store a reference to a matrix | ||
38 | 2857 | explicit ParallelMultiTypeMatrixAdapter (const M& A) | |
39 |
2/4✓ Branch 1 taken 2844 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2844 times.
✗ Branch 5 not taken.
|
2857 | : ParentType(A) {} |
40 | |||
41 | //! constructor: store an std::shared_ptr to a matrix | ||
42 | explicit ParallelMultiTypeMatrixAdapter (std::shared_ptr<const M> A) | ||
43 | : ParentType(A) {} | ||
44 | |||
45 | //! apply operator to x: \f$ y = A(x) \f$ | ||
46 | 162742 | void apply (const X& x, Y& y) const override | |
47 | { | ||
48 | 162742 | y = 0.0; | |
49 | |||
50 | 162742 | auto& A = this->getmat(); | |
51 | using namespace Dune::Hybrid; | ||
52 | 325484 | forEach(integralRange(Dune::Hybrid::size(y)), [&](auto&& i) | |
53 | { | ||
54 | 162742 | forEach(integralRange(Dune::Hybrid::size(x)), [&](auto&& j) | |
55 | { | ||
56 | 325484 | const auto& mat = A[i][j]; | |
57 | 650968 | const auto& xx = x[j]; | |
58 | 325484 | auto& yy = y[i]; | |
59 | |||
60 | 8051364286 | Dumux::parallelFor(mat.N(), [&](const std::size_t ii) | |
61 | { | ||
62 | 3167918110 | const auto& row = mat[ii]; | |
63 | 3167918110 | const auto endj = row.end(); | |
64 |
8/8✓ Branch 0 taken 644188900 times.
✓ Branch 1 taken 202538082 times.
✓ Branch 2 taken 391374018 times.
✓ Branch 3 taken 202487683 times.
✓ Branch 4 taken 724172037 times.
✓ Branch 5 taken 1381171473 times.
✓ Branch 6 taken 2901225649 times.
✓ Branch 7 taken 1381720872 times.
|
4334754565 | for (auto jj=row.begin(); jj!=endj; ++jj) |
65 | { | ||
66 | 3477022755 | const auto& xj = Dune::Impl::asVector(xx[jj.index()]); | |
67 | 3477022755 | auto&& yi = Dune::Impl::asVector(yy[ii]); | |
68 | 6954045510 | Dune::Impl::asMatrix(*jj).umv(xj, yi); | |
69 | } | ||
70 | }); | ||
71 | }); | ||
72 | }); | ||
73 | 162742 | } | |
74 | |||
75 | //! apply operator to x, scale and add: \f$ y = y + \alpha A(x) \f$ | ||
76 | 2956 | void applyscaleadd (field_type alpha, const X& x, Y& y) const override | |
77 | { | ||
78 | 2956 | auto& A = this->getmat(); | |
79 | using namespace Dune::Hybrid; | ||
80 | 26604 | forEach(integralRange(Dune::Hybrid::size(y)), [&](auto&& i) | |
81 | { | ||
82 | 5912 | forEach(integralRange(Dune::Hybrid::size(x)), [&](auto&& j) | |
83 | { | ||
84 | 5912 | const auto& mat = A[i][j]; | |
85 | 11824 | const auto& xx = x[j]; | |
86 | 5912 | auto& yy = y[i]; | |
87 | |||
88 | 122210392 | Dumux::parallelFor(mat.N(), [&](const std::size_t ii) | |
89 | { | ||
90 | 30551120 | const auto& row = mat[ii]; | |
91 | 30551120 | const auto endj = row.end(); | |
92 |
8/8✓ Branch 0 taken 5822692 times.
✓ Branch 1 taken 1902922 times.
✓ Branch 2 taken 3199214 times.
✓ Branch 3 taken 1852422 times.
✓ Branch 4 taken 6113756 times.
✓ Branch 5 taken 13122638 times.
✓ Branch 6 taken 80191180 times.
✓ Branch 7 taken 13673138 times.
|
125877962 | for (auto jj=row.begin(); jj!=endj; ++jj) |
93 | { | ||
94 | 95326842 | const auto& xj = Dune::Impl::asVector(xx[jj.index()]); | |
95 | 95326842 | auto&& yi = Dune::Impl::asVector(yy[ii]); | |
96 | 190653684 | Dune::Impl::asMatrix(*jj).usmv(alpha, xj, yi); | |
97 | } | ||
98 | }); | ||
99 | }); | ||
100 | }); | ||
101 | 2956 | } | |
102 | }; | ||
103 | |||
104 | } // end namespace Dumux | ||
105 | |||
106 | #endif | ||
107 |