21#ifndef OPM_PARALLELISTLINFORMATION_HEADER_INCLUDED
22#define OPM_PARALLELISTLINFORMATION_HEADER_INCLUDED
26#if HAVE_MPI && HAVE_DUNE_ISTL
34#include <dune/istl/owneroverlapcopy.hh>
43class ParallelISTLInformation
47 using ParallelIndexSet = Dune::OwnerOverlapCopyCommunication<int, int>::ParallelIndexSet;
49 using RemoteIndices = Dune::OwnerOverlapCopyCommunication<int, int>::RemoteIndices;
52 ParallelISTLInformation();
56 explicit ParallelISTLInformation(MPI_Comm communicator);
62 ParallelISTLInformation(
const std::shared_ptr<ParallelIndexSet>& indexSet,
63 const std::shared_ptr<RemoteIndices>& remoteIndices,
64 MPI_Comm communicator);
69 ParallelISTLInformation(
const ParallelISTLInformation& other);
72 std::shared_ptr<ParallelIndexSet> indexSet()
const
78 std::shared_ptr<RemoteIndices> remoteIndices()
const
80 return remoteIndices_;
91 void copyValuesTo(ParallelIndexSet& indexSet, RemoteIndices& remoteIndices,
92 std::size_t local_component_size = 0,
93 std::size_t num_components = 1)
const;
99 void copyOwnerToAll (
const T& source, T& dest)
const;
102 const std::vector<double>& updateOwnerMask(
const T& container)
const;
109 const std::vector<double>& getOwnerMask()
const
133 template<
typename Container,
typename BinaryOperator,
typename T>
134 void computeReduction(
const Container& container, BinaryOperator binaryOperator,
138 template<
typename... Containers,
typename... BinaryOperators,
typename... ReturnValues>
139 void computeTupleReduction(
const std::tuple<Containers...>& containers,
140 std::tuple<BinaryOperators...>& operators,
141 std::tuple<ReturnValues...>& values)
const;
143 std::shared_ptr<ParallelIndexSet> indexSet_;
144 std::shared_ptr<RemoteIndices> remoteIndices_;
146 mutable std::vector<double> ownerMask_;
156 template<
typename BinaryOperator>
157 struct MaskIDOperator
165 template<
class T,
class T1>
166 T operator()(
const T& t1,
const T& t2,
const T1& mask)
168 return b_(t1, maskValue(t2, mask));
170 template<
class T,
class T1>
171 T maskValue(
const T& t,
const T1& mask)
175 BinaryOperator& localOperator()
190 struct InnerProductFunctor
199 T operator()(
const T& t1,
const T& t2,
const T1& mask)
201 T masked = maskValue(t2, mask);
202 return t1 + masked * masked;
205 T maskValue(
const T& t,
const T1& mask)
209 std::plus<T> localOperator()
211 return std::plus<T>();
213 template<
class V, std::enable_if_t<std::is_convertible_v<T, V>,
int> = 0>
225 template<
typename BinaryOperator>
226 struct MaskToMinOperator
228 explicit MaskToMinOperator(BinaryOperator b)
237 template<
class T,
class T1>
238 T operator()(
const T& t1,
const T& t2,
const T1& mask)
240 return b_(t1, maskValue(t2, mask));
242 template<
class T,
class T1>
243 T maskValue(
const T& t,
const T1& mask)
251 return getInitialValue<T>();
257 return std::numeric_limits<T>::lowest();
263 BinaryOperator& localOperator()
274 template<
typename BinaryOperator>
275 struct MaskToMaxOperator
277 explicit MaskToMaxOperator(BinaryOperator b)
286 template<
class T,
class T1>
287 T operator()(
const T& t1,
const T& t2,
const T1& mask)
289 return b_(t1, maskValue(t2, mask));
291 template<
class T,
class T1>
292 T maskValue(
const T& t,
const T1& mask)
300 return std::numeric_limits<T>::max();
303 BinaryOperator& localOperator()
310 return std::numeric_limits<T>::max();
319 MaskIDOperator<std::plus<T> >
320 makeGlobalSumFunctor()
322 return MaskIDOperator<std::plus<T> >();
328 auto makeGlobalMaxFunctor()
332 using result_type = T;
333 const result_type& operator()(
const T& t1,
const T& t2)
335 return std::max(t1, t2);
338 return MaskToMinOperator(MaxOp());
344 template<
typename T,
typename Enable =
void>
347 using result_type = T;
348 result_type operator()(
const T& t1,
351 return std::max(std::abs(t1), std::abs(t2));
359 struct MaxAbsFunctor<T, typename std::enable_if<std::is_unsigned<T>::value>::type>
361 using result_type = T;
362 result_type operator()(
const T& t1,
365 return std::max(t1, t2);
374 MaskIDOperator<detail::MaxAbsFunctor<T> >
375 makeLInfinityNormFunctor()
377 return MaskIDOperator<detail::MaxAbsFunctor<T> >();
384 makeGlobalMinFunctor()
388 using result_type = T;
389 const result_type& operator()(
const T& t1,
const T& t2)
391 return std::min(t1, t2);
394 return MaskToMaxOperator(MinOp());
397 InnerProductFunctor<T>
398 makeInnerProductFunctor()
400 return InnerProductFunctor<T>();
418 ->
decltype(container[0]*(*maskContainer)[0]);
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilbioeffectsmodules.hh:43
auto accumulateMaskedValues(const T1 &container, const std::vector< double > *maskContainer) -> decltype(container[0] *(*maskContainer)[0])
Accumulates entries masked with 1.