ParallelIstlInformation.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014, 2015 Dr. Markus Blatt - HPC-Simulation-Software & Services
3  Copyright 2014, 2015 Statoil ASA
4  Copyright 2015 NTNU
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 3 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
21 #ifndef OPM_PARALLELISTLINFORMTION_HEADER_INCLUDED
22 #define OPM_PARALLELISTLINFORMTION_HEADER_INCLUDED
23 
24 
25 #include <opm/core/grid.h>
26 #include <opm/common/ErrorMacros.hpp>
27 #include <boost/any.hpp>
28 #include <exception>
29 
30 #if HAVE_MPI && HAVE_DUNE_ISTL
31 
32 #include <opm/common/utility/platform_dependent/disable_warnings.h>
33 #include <mpi.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>
39 
40 #include <algorithm>
41 #include <functional>
42 #include <limits>
43 #include <numeric>
44 #include <type_traits>
45 
46 namespace Opm
47 {
48 namespace
49 {
50 
51  template<class T>
52  struct is_tuple
53  : std::integral_constant<bool, false>
54  {};
55  template<typename... T>
56  struct is_tuple<std::tuple<T...> >
57  : std::integral_constant<bool, true>
58  {};
59 }
60 
63 class ParallelISTLInformation
64 {
65 public:
67  typedef Dune::OwnerOverlapCopyCommunication<int, int>::ParallelIndexSet ParallelIndexSet;
69  typedef Dune::OwnerOverlapCopyCommunication<int, int>::RemoteIndices RemoteIndices;
70 
72  ParallelISTLInformation()
73  : indexSet_(new ParallelIndexSet),
74  remoteIndices_(new RemoteIndices(*indexSet_, *indexSet_, MPI_COMM_WORLD)),
75  communicator_(MPI_COMM_WORLD)
76  {}
79  ParallelISTLInformation(MPI_Comm communicator)
80  : indexSet_(new ParallelIndexSet),
81  remoteIndices_(new RemoteIndices(*indexSet_, *indexSet_, communicator)),
82  communicator_(communicator)
83  {}
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)
92  {}
96  ParallelISTLInformation(const ParallelISTLInformation& other)
97  : indexSet_(other.indexSet_), remoteIndices_(other.remoteIndices_),
98  communicator_(other.communicator_)
99  {}
101  std::shared_ptr<ParallelIndexSet> indexSet() const
102  {
103  return indexSet_;
104  }
106  std::shared_ptr<RemoteIndices> remoteIndices() const
107  {
108  return remoteIndices_;
109  }
111  Dune::CollectiveCommunication<MPI_Comm> communicator() const
112  {
113  return communicator_;
114  }
118  void copyValuesTo(ParallelIndexSet& indexSet, RemoteIndices& remoteIndices,
119  std::size_t local_component_size = 0, std::size_t num_components = 1) const
120  {
121  ParallelIndexSet::GlobalIndex global_component_size = local_component_size;
122  if ( num_components > 1 )
123  {
124  ParallelIndexSet::GlobalIndex max_gi = 0;
125  // component the max global index
126  for( auto i = indexSet_->begin(), end = indexSet_->end(); i != end; ++i )
127  {
128  max_gi = std::max(max_gi, i->global());
129  }
130  global_component_size = max_gi+1;
131  global_component_size = communicator_.max(global_component_size);
132  }
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>();
139  }
143  template<class T>
144  void copyOwnerToAll (const T& source, T& dest) const
145  {
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;
150  AllSet destFlags;
151  Dune::Interface interface(communicator_);
152  if( !remoteIndices_->isSynced() )
153  {
154  remoteIndices_->rebuild<false>();
155  }
156  interface.build(*remoteIndices_,sourceFlags,destFlags);
157  Dune::BufferedCommunicator communicator;
158  communicator.template build<T>(interface);
159  communicator.template forward<CopyGatherScatter<T> >(source,dest);
160  communicator.free();
161  }
162  template<class T>
163  void updateOwnerMask(const T& container) const
164  {
165  if( ! indexSet_ )
166  {
167  OPM_THROW(std::runtime_error, "Trying to update owner mask without parallel information!");
168  }
169  if( static_cast<std::size_t>(container.size())!= ownerMask_.size() )
170  {
171  ownerMask_.resize(container.size(), 1.);
172  for( auto i=indexSet_->begin(), end=indexSet_->end(); i!=end; ++i )
173  {
174  if (i->local().attribute()!=Dune::OwnerOverlapCopyAttributeSet::owner)
175  {
176  ownerMask_[i->local().local()] = 0.;
177  }
178  }
179  }
180  }
199  template<typename Container, typename BinaryOperator, typename T>
200  void computeReduction(const Container& container, BinaryOperator binaryOperator,
201  T& value) const
202  {
203  computeReduction(container, binaryOperator, value, is_tuple<Container>());
204  }
205 private:
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
212  {
213  computeTupleReduction(container, binaryOperator, value);
214  }
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
221  {
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);
227  }
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
233  {
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 )
241  {
242  return;
243  }
244  // Copy the initial values.
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]));
250  values=init;
251  for( auto rvals=receivedValues.begin(), endvals=receivedValues.end(); rvals!=endvals;
252  ++rvals )
253  {
254  computeGlobalReduction(*rvals, operators, values);
255  }
256  }
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
265  {}
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
272  {
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);
276  }
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
285  {}
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
292  {
293  const auto& container = std::get<I>(containers);
294  if( container.size() )
295  {
296  auto& reduceOperator = std::get<I>(operators);
297  // Eigen:Block does not support STL iterators!!!!
298  // Therefore we need to rely on the harder random-access
299  // property of the containers. But this should be save, too.
300  // Just commenting out code in the hope that Eigen might improve
301  // in this regard in the future.
302  //auto newVal = container.begin();
303  auto mask = ownerMask_.begin();
304  auto& value = std::get<I>(values);
305  value = reduceOperator.getInitialValue();
306 
307  for( auto endVal=ownerMask_.end(); mask!=endVal;
308  /*++newVal,*/ ++mask )
309  {
310  value = reduceOperator(value, container[mask-ownerMask_.begin()], *mask);
311  }
312  }
313  computeLocalReduction<I+1>(containers, operators, values);
314  }
316  template<typename T>
317  struct CopyGatherScatter
318  {
319  typedef typename Dune::CommPolicy<T>::IndexedType V;
320 
321  static V gather(const T& a, std::size_t i)
322  {
323  return a[i];
324  }
325 
326  static void scatter(T& a, V v, std::size_t i)
327  {
328  a[i] = v;
329  }
330  };
331  template<class T>
332  class IndexSetInserter
333  {
334  public:
335  typedef T ParallelIndexSet;
336  typedef typename ParallelIndexSet::LocalIndex LocalIndex;
337  typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;
338 
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)
344  {}
345  void operator()(const typename ParallelIndexSet::IndexPair& pair)
346  {
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()));
351  }
352  private:
353  ParallelIndexSet* indexSet_;
355  GlobalIndex component_size_;
357  std::size_t local_component_size_;
359  std::size_t num_components_;
360  };
361  std::shared_ptr<ParallelIndexSet> indexSet_;
362  std::shared_ptr<RemoteIndices> remoteIndices_;
363  Dune::CollectiveCommunication<MPI_Comm> communicator_;
364  mutable std::vector<double> ownerMask_;
365 };
366 
367  namespace Reduction
368  {
373  // the reduction operation.
374  template<typename BinaryOperator>
375  struct MaskIDOperator
376  {
377  // This is a real nice one: numeric limits needs a type without const
378  // or reference qualifier. Otherwise we get complete nonesense.
379  typedef typename std::remove_cv<
380  typename std::remove_reference<typename BinaryOperator::result_type>::type
381  >::type Result;
388  template<class T, class T1>
389  T operator()(const T& t1, const T& t2, const T1& mask)
390  {
391  return b_(t1, maskValue(t2, mask));
392  }
393  template<class T, class T1>
394  T maskValue(const T& t, const T1& mask)
395  {
396  return t*mask;
397  }
398  BinaryOperator& localOperator()
399  {
400  return b_;
401  }
402  Result getInitialValue()
403  {
404  return Result();
405  }
406  private:
407  BinaryOperator b_;
408  };
409 
411  template<class T>
412  struct InnerProductFunctor
413  {
420  template<class T1>
421  T operator()(const T& t1, const T& t2, const T1& mask)
422  {
423  T masked = maskValue(t2, mask);
424  return t1 + masked * masked;
425  }
426  template<class T1>
427  T maskValue(const T& t, const T1& mask)
428  {
429  return t*mask;
430  }
431  std::plus<T> localOperator()
432  {
433  return std::plus<T>();
434  }
435  T getInitialValue()
436  {
437  return T();
438  }
439  };
440 
445  // the reduction operation.
446  template<typename BinaryOperator>
447  struct MaskToMinOperator
448  {
449  // This is a real nice one: numeric limits has to a type without const
450  // or reference. Otherwise we get complete nonesense.
451  typedef typename std::remove_reference<
452  typename std::remove_const<typename BinaryOperator::result_type>::type
453  >::type Result;
454 
455  MaskToMinOperator(BinaryOperator b)
456  : b_(b)
457  {}
464  template<class T, class T1>
465  T operator()(const T& t1, const T& t2, const T1& mask)
466  {
467  return b_(t1, maskValue(t2, mask));
468  }
469  template<class T, class T1>
470  T maskValue(const T& t, const T1& mask)
471  {
472  if( mask )
473  {
474  return t;
475  }
476  else
477  {
478  return getInitialValue();
479  }
480  }
481  Result getInitialValue()
482  {
483  //g++-4.4 does not support std::numeric_limits<T>::lowest();
484  // we rely on IEE 754 for floating point values and use min()
485  // for integral types.
486  if( std::is_integral<Result>::value )
487  {
488  return std::numeric_limits<Result>::min();
489  }
490  else
491  {
492  return -std::numeric_limits<Result>::max();
493  }
494  }
499  BinaryOperator& localOperator()
500  {
501  return b_;
502  }
503  private:
504  BinaryOperator b_;
505  };
506 
510  template<typename BinaryOperator>
511  struct MaskToMaxOperator
512  {
513  // This is a real nice one: numeric limits has to a type without const
514  // or reference. Otherwise we get complete nonesense.
515  typedef typename std::remove_cv<
516  typename std::remove_reference<typename BinaryOperator::result_type>::type
517  >::type Result;
518 
519  MaskToMaxOperator(BinaryOperator b)
520  : b_(b)
521  {}
528  template<class T, class T1>
529  T operator()(const T& t1, const T& t2, const T1& mask)
530  {
531  return b_(t1, maskValue(t2, mask));
532  }
533  template<class T, class T1>
534  T maskValue(const T& t, const T1& mask)
535  {
536  if( mask )
537  {
538  return t;
539  }
540  else
541  {
542  return std::numeric_limits<T>::max();
543  }
544  }
545  BinaryOperator& localOperator()
546  {
547  return b_;
548  }
549  Result getInitialValue()
550  {
551  return std::numeric_limits<Result>::max();
552  }
553  private:
554  BinaryOperator b_;
555  };
559  template<class T>
560  MaskIDOperator<std::plus<T> >
561  makeGlobalSumFunctor()
562  {
563  return MaskIDOperator<std::plus<T> >();
564  }
568  template<class T>
569  MaskToMinOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
570  makeGlobalMaxFunctor()
571  {
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>));
575  }
579  template<class T>
580  MaskToMaxOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
581  makeGlobalMinFunctor()
582  {
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>));
586  }
587  template<class T>
588  InnerProductFunctor<T>
589  makeInnerProductFunctor()
590  {
591  return InnerProductFunctor<T>();
592  }
593  } // end namespace Reduction
594 } // end namespace Opm
595 
596 #endif
597 
598 namespace Opm
599 {
609 
610 inline void extractParallelGridInformationToISTL(boost::any& anyComm, const UnstructuredGrid& grid)
611 {
612  (void)anyComm; (void)grid;
613 }
614 } // end namespace Opm
615 
616 #endif
Definition: grid.h:98
Definition: AnisotropicEikonal.hpp:43
STL namespace.
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