12#ifndef DUMUX_LINEAR_PARALLEL_MATRIX_ADAPTER_HH
13#define DUMUX_LINEAR_PARALLEL_MATRIX_ADAPTER_HH
15#include <dune/common/hybridutilities.hh>
16#include <dune/istl/operators.hh>
17#include <dune/istl/bcrsmatrix.hh>
26template<
class M,
class X,
class Y>
29 using ParentType = Dune::MatrixAdapter<M,X,Y>;
46 void apply (
const X& x, Y& y)
const override
50 auto& A = this->getmat();
51 using namespace Dune::Hybrid;
52 forEach(integralRange(Dune::Hybrid::size(y)), [&](
auto&& i)
54 forEach(integralRange(Dune::Hybrid::size(x)), [&](
auto&& j)
56 const auto& mat = A[i][j];
57 const auto& xx = x[j];
62 const auto& row = mat[ii];
63 const auto endj = row.end();
64 for (auto jj=row.begin(); jj!=endj; ++jj)
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);
78 auto& A = this->getmat();
79 using namespace Dune::Hybrid;
80 forEach(integralRange(Dune::Hybrid::size(y)), [&](
auto&& i)
82 forEach(integralRange(Dune::Hybrid::size(x)), [&](
auto&& j)
84 const auto& mat = A[i][j];
85 const auto& xx = x[j];
90 const auto& row = mat[ii];
91 const auto endj = row.end();
92 for (auto jj=row.begin(); jj!=endj; ++jj)
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);
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
Parallel for loop (multithreading)