21 #ifndef OPM_PARALLELISTLINFORMTION_HEADER_INCLUDED
22 #define OPM_PARALLELISTLINFORMTION_HEADER_INCLUDED
26 #include <opm/common/ErrorMacros.hpp>
27 #include <boost/any.hpp>
30 #if HAVE_MPI && HAVE_DUNE_ISTL
32 #include <opm/common/utility/platform_dependent/disable_warnings.h>
34 #include <dune/istl/owneroverlapcopy.hh>
35 #include <dune/common/parallel/interface.hh>
36 #include <dune/common/parallel/communicator.hh>
37 #include <dune/common/enumset.hh>
38 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
44 #include <type_traits>
53 : std::integral_constant<bool, false>
55 template<
typename... T>
56 struct is_tuple<
std::tuple<T...> >
57 : std::integral_constant<bool, true>
63 class ParallelISTLInformation
67 typedef Dune::OwnerOverlapCopyCommunication<int, int>::ParallelIndexSet ParallelIndexSet;
69 typedef Dune::OwnerOverlapCopyCommunication<int, int>::RemoteIndices RemoteIndices;
72 ParallelISTLInformation()
73 : indexSet_(new ParallelIndexSet),
74 remoteIndices_(new RemoteIndices(*indexSet_, *indexSet_, MPI_COMM_WORLD)),
75 communicator_(MPI_COMM_WORLD)
79 ParallelISTLInformation(MPI_Comm communicator)
80 : indexSet_(new ParallelIndexSet),
81 remoteIndices_(new RemoteIndices(*indexSet_, *indexSet_, communicator)),
82 communicator_(communicator)
88 ParallelISTLInformation(
const std::shared_ptr<ParallelIndexSet>& indexSet,
89 const std::shared_ptr<RemoteIndices>& remoteIndices,
90 MPI_Comm communicator)
91 : indexSet_(indexSet), remoteIndices_(remoteIndices), communicator_(communicator)
96 ParallelISTLInformation(
const ParallelISTLInformation& other)
97 : indexSet_(other.indexSet_), remoteIndices_(other.remoteIndices_),
98 communicator_(other.communicator_)
101 std::shared_ptr<ParallelIndexSet> indexSet()
const
106 std::shared_ptr<RemoteIndices> remoteIndices()
const
108 return remoteIndices_;
111 Dune::CollectiveCommunication<MPI_Comm> communicator()
const
113 return communicator_;
118 void copyValuesTo(ParallelIndexSet& indexSet, RemoteIndices& remoteIndices,
119 std::size_t local_component_size = 0, std::size_t num_components = 1)
const
121 ParallelIndexSet::GlobalIndex global_component_size = local_component_size;
122 if ( num_components > 1 )
124 ParallelIndexSet::GlobalIndex max_gi = 0;
126 for(
auto i = indexSet_->begin(), end = indexSet_->end(); i != end; ++i )
128 max_gi = std::max(max_gi, i->global());
130 global_component_size = max_gi+1;
131 global_component_size = communicator_.max(global_component_size);
133 indexSet.beginResize();
134 IndexSetInserter<ParallelIndexSet> inserter(indexSet, global_component_size,
135 local_component_size, num_components);
136 std::for_each(indexSet_->begin(), indexSet_->end(), inserter);
137 indexSet.endResize();
138 remoteIndices.rebuild<
false>();
144 void copyOwnerToAll (
const T& source, T& dest)
const
146 typedef Dune::Combine<Dune::EnumItem<Dune::OwnerOverlapCopyAttributeSet::AttributeSet,Dune::OwnerOverlapCopyAttributeSet::owner>,Dune::EnumItem<Dune::OwnerOverlapCopyAttributeSet::AttributeSet,Dune::OwnerOverlapCopyAttributeSet::overlap>,Dune::OwnerOverlapCopyAttributeSet::AttributeSet> OwnerOverlapSet;
147 typedef Dune::EnumItem<Dune::OwnerOverlapCopyAttributeSet::AttributeSet,Dune::OwnerOverlapCopyAttributeSet::owner> OwnerSet;
148 typedef Dune::Combine<OwnerOverlapSet, Dune::EnumItem<Dune::OwnerOverlapCopyAttributeSet::AttributeSet,Dune::OwnerOverlapCopyAttributeSet::copy>,Dune::OwnerOverlapCopyAttributeSet::AttributeSet> AllSet;
149 OwnerSet sourceFlags;
151 Dune::Interface interface(communicator_);
152 if( !remoteIndices_->isSynced() )
154 remoteIndices_->rebuild<
false>();
156 interface.build(*remoteIndices_,sourceFlags,destFlags);
157 Dune::BufferedCommunicator communicator;
158 communicator.template build<T>(interface);
159 communicator.template forward<CopyGatherScatter<T> >(source,dest);
163 void updateOwnerMask(
const T& container)
const
167 OPM_THROW(std::runtime_error,
"Trying to update owner mask without parallel information!");
169 if( static_cast<std::size_t>(container.size())!= ownerMask_.size() )
171 ownerMask_.resize(container.size(), 1.);
172 for(
auto i=indexSet_->begin(), end=indexSet_->end(); i!=end; ++i )
174 if (i->local().attribute()!=Dune::OwnerOverlapCopyAttributeSet::owner)
176 ownerMask_[i->local().local()] = 0.;
199 template<
typename Container,
typename BinaryOperator,
typename T>
200 void computeReduction(
const Container& container, BinaryOperator binaryOperator,
203 computeReduction(container, binaryOperator, value, is_tuple<Container>());
209 template<
typename Container,
typename BinaryOperator,
typename T>
210 void computeReduction(
const Container& container, BinaryOperator binaryOperator,
211 T& value, std::integral_constant<bool,true>)
const
213 computeTupleReduction(container, binaryOperator, value);
218 template<
typename Container,
typename BinaryOperator,
typename T>
219 void computeReduction(
const Container& container, BinaryOperator binaryOperator,
220 T& value, std::integral_constant<bool,false>)
const
222 std::tuple<const Container&> containers=std::tuple<const Container&>(container);
223 auto values=std::make_tuple(value);
224 auto operators=std::make_tuple(binaryOperator);
225 computeTupleReduction(containers, operators, values);
226 value=std::get<0>(values);
229 template<
typename... Containers,
typename... BinaryOperators,
typename... ReturnValues>
230 void computeTupleReduction(
const std::tuple<Containers...>& containers,
231 std::tuple<BinaryOperators...>& operators,
232 std::tuple<ReturnValues...>& values)
const
234 static_assert(std::tuple_size<std::tuple<Containers...> >::value==
235 std::tuple_size<std::tuple<BinaryOperators...> >::value,
236 "We need the same number of containers and binary operators");
237 static_assert(std::tuple_size<std::tuple<Containers...> >::value==
238 std::tuple_size<std::tuple<ReturnValues...> >::value,
239 "We need the same number of containers and return values");
240 if( std::tuple_size<std::tuple<Containers...> >::value==0 )
245 std::tuple<ReturnValues...> init=values;
246 updateOwnerMask(std::get<0>(containers));
247 computeLocalReduction(containers, operators, values);
248 std::vector<std::tuple<ReturnValues...> > receivedValues(communicator_.size());
249 communicator_.allgather(&values, 1, &(receivedValues[0]));
251 for(
auto rvals=receivedValues.begin(), endvals=receivedValues.end(); rvals!=endvals;
254 computeGlobalReduction(*rvals, operators, values);
260 template<
int I=0,
typename... BinaryOperators,
typename... ReturnValues>
261 typename std::enable_if<I ==
sizeof...(BinaryOperators),
void>::type
262 computeGlobalReduction(
const std::tuple<ReturnValues...>&,
263 std::tuple<BinaryOperators...>&,
264 std::tuple<ReturnValues...>&)
const
267 template<
int I=0,
typename... BinaryOperators,
typename... ReturnValues>
268 typename std::enable_if<I !=
sizeof...(BinaryOperators),
void>::type
269 computeGlobalReduction(
const std::tuple<ReturnValues...>& receivedValues,
270 std::tuple<BinaryOperators...>& operators,
271 std::tuple<ReturnValues...>& values)
const
273 auto& val=std::get<I>(values);
274 val = std::get<I>(operators).localOperator()(val, std::get<I>(receivedValues));
275 computeGlobalReduction<I+1>(receivedValues, operators, values);
280 template<
int I=0,
typename... Containers,
typename... BinaryOperators,
typename... ReturnValues>
281 typename std::enable_if<I==
sizeof...(Containers),
void>::type
282 computeLocalReduction(
const std::tuple<Containers...>&,
283 std::tuple<BinaryOperators...>&,
284 std::tuple<ReturnValues...>&)
const
287 template<
int I=0,
typename... Containers,
typename... BinaryOperators,
typename... ReturnValues>
288 typename std::enable_if<I!=
sizeof...(Containers),
void>::type
289 computeLocalReduction(
const std::tuple<Containers...>& containers,
290 std::tuple<BinaryOperators...>& operators,
291 std::tuple<ReturnValues...>& values)
const
293 const auto& container = std::get<I>(containers);
294 if( container.size() )
296 auto& reduceOperator = std::get<I>(operators);
303 auto mask = ownerMask_.begin();
304 auto& value = std::get<I>(values);
305 value = reduceOperator.getInitialValue();
307 for(
auto endVal=ownerMask_.end(); mask!=endVal;
310 value = reduceOperator(value, container[mask-ownerMask_.begin()], *mask);
313 computeLocalReduction<I+1>(containers, operators, values);
317 struct CopyGatherScatter
319 typedef typename Dune::CommPolicy<T>::IndexedType V;
321 static V gather(
const T& a, std::size_t i)
326 static void scatter(T& a, V v, std::size_t i)
332 class IndexSetInserter
335 typedef T ParallelIndexSet;
336 typedef typename ParallelIndexSet::LocalIndex LocalIndex;
337 typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
339 IndexSetInserter(ParallelIndexSet& indexSet,
const GlobalIndex& component_size,
340 std::size_t local_component_size, std::size_t num_components)
341 : indexSet_(&indexSet), component_size_(component_size),
342 local_component_size_(local_component_size),
343 num_components_(num_components)
345 void operator()(
const typename ParallelIndexSet::IndexPair& pair)
347 for(std::size_t i = 0; i < num_components_; i++)
348 indexSet_->add(i * component_size_ + pair.global(),
349 LocalIndex(i * local_component_size_ + pair.local(),
350 pair.local().attribute()));
353 ParallelIndexSet* indexSet_;
355 GlobalIndex component_size_;
357 std::size_t local_component_size_;
359 std::size_t num_components_;
361 std::shared_ptr<ParallelIndexSet> indexSet_;
362 std::shared_ptr<RemoteIndices> remoteIndices_;
363 Dune::CollectiveCommunication<MPI_Comm> communicator_;
364 mutable std::vector<double> ownerMask_;
374 template<
typename BinaryOperator>
375 struct MaskIDOperator
379 typedef typename std::remove_cv<
380 typename std::remove_reference<typename BinaryOperator::result_type>::type
388 template<
class T,
class T1>
389 T operator()(
const T& t1,
const T& t2,
const T1& mask)
391 return b_(t1, maskValue(t2, mask));
393 template<
class T,
class T1>
394 T maskValue(
const T& t,
const T1& mask)
398 BinaryOperator& localOperator()
402 Result getInitialValue()
412 struct InnerProductFunctor
421 T operator()(
const T& t1,
const T& t2,
const T1& mask)
423 T masked = maskValue(t2, mask);
424 return t1 + masked * masked;
427 T maskValue(
const T& t,
const T1& mask)
431 std::plus<T> localOperator()
433 return std::plus<T>();
446 template<
typename BinaryOperator>
447 struct MaskToMinOperator
451 typedef typename std::remove_reference<
452 typename std::remove_const<typename BinaryOperator::result_type>::type
455 MaskToMinOperator(BinaryOperator b)
464 template<
class T,
class T1>
465 T operator()(
const T& t1,
const T& t2,
const T1& mask)
467 return b_(t1, maskValue(t2, mask));
469 template<
class T,
class T1>
470 T maskValue(
const T& t,
const T1& mask)
478 return getInitialValue();
481 Result getInitialValue()
486 if( std::is_integral<Result>::value )
488 return std::numeric_limits<Result>::min();
492 return -std::numeric_limits<Result>::max();
499 BinaryOperator& localOperator()
510 template<
typename BinaryOperator>
511 struct MaskToMaxOperator
515 typedef typename std::remove_cv<
516 typename std::remove_reference<typename BinaryOperator::result_type>::type
519 MaskToMaxOperator(BinaryOperator b)
528 template<
class T,
class T1>
529 T operator()(
const T& t1,
const T& t2,
const T1& mask)
531 return b_(t1, maskValue(t2, mask));
533 template<
class T,
class T1>
534 T maskValue(
const T& t,
const T1& mask)
542 return std::numeric_limits<T>::max();
545 BinaryOperator& localOperator()
549 Result getInitialValue()
551 return std::numeric_limits<Result>::max();
560 MaskIDOperator<std::plus<T> >
561 makeGlobalSumFunctor()
563 return MaskIDOperator<std::plus<T> >();
569 MaskToMinOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
570 makeGlobalMaxFunctor()
572 return MaskToMinOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
573 (std::pointer_to_binary_function<const T&,const T&,const T&>
574 ((
const T&(*)(
const T&,
const T&))
std::max<T>));
580 MaskToMaxOperator<
std::pointer_to_binary_function<const T&,const T&,const T&> >
581 makeGlobalMinFunctor()
583 return MaskToMaxOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
584 (std::pointer_to_binary_function<const T&,const T&,const T&>
585 ((
const T&(*)(
const T&,
const T&))
std::min<T>));
588 InnerProductFunctor<T>
589 makeInnerProductFunctor()
591 return InnerProductFunctor<T>();
612 (void)anyComm; (void)grid;
Definition: AnisotropicEikonal.hpp:43
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
void extractParallelGridInformationToISTL(boost::any &anyComm, const UnstructuredGrid &grid)
Extracts the information about the data decomposition from the grid for dune-istl.
Definition: ParallelIstlInformation.hpp:610