Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
partitioner.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
18 #include <deal.II/base/partitioner.templates.h>
19 
21 
22 namespace Utilities
23 {
24  namespace MPI
25  {
27  : global_size(0)
28  , local_range_data(
29  std::pair<types::global_dof_index, types::global_dof_index>(0, 0))
30  , n_ghost_indices_data(0)
31  , n_import_indices_data(0)
32  , n_ghost_indices_in_larger_set(0)
33  , my_pid(0)
34  , n_procs(1)
35  , communicator(MPI_COMM_SELF)
36  , have_ghost_indices(false)
37  {}
38 
39 
40 
41  Partitioner::Partitioner(const unsigned int size)
42  : global_size(size)
43  , locally_owned_range_data(size)
44  , local_range_data(
45  std::pair<types::global_dof_index, types::global_dof_index>(0, size))
46  , n_ghost_indices_data(0)
47  , n_import_indices_data(0)
48  , n_ghost_indices_in_larger_set(0)
49  , my_pid(0)
50  , n_procs(1)
51  , communicator(MPI_COMM_SELF)
52  , have_ghost_indices(false)
53  {
57  }
58 
59 
60 
61  Partitioner::Partitioner(const IndexSet &locally_owned_indices,
62  const IndexSet &ghost_indices_in,
63  const MPI_Comm communicator_in)
64  : global_size(
65  static_cast<types::global_dof_index>(locally_owned_indices.size()))
66  , n_ghost_indices_data(0)
67  , n_import_indices_data(0)
68  , n_ghost_indices_in_larger_set(0)
69  , my_pid(0)
70  , n_procs(1)
71  , communicator(communicator_in)
72  , have_ghost_indices(false)
73  {
74  set_owned_indices(locally_owned_indices);
75  set_ghost_indices(ghost_indices_in);
76  }
77 
78 
79 
80  Partitioner::Partitioner(const IndexSet &locally_owned_indices,
81  const MPI_Comm communicator_in)
82  : global_size(
83  static_cast<types::global_dof_index>(locally_owned_indices.size()))
84  , n_ghost_indices_data(0)
85  , n_import_indices_data(0)
86  , n_ghost_indices_in_larger_set(0)
87  , my_pid(0)
88  , n_procs(1)
89  , communicator(communicator_in)
90  , have_ghost_indices(false)
91  {
92  set_owned_indices(locally_owned_indices);
93  }
94 
95 
96 
97  void
98  Partitioner::reinit(const IndexSet &vector_space_vector_index_set,
99  const IndexSet &read_write_vector_index_set,
100  const MPI_Comm &communicator_in)
101  {
102  have_ghost_indices = false;
103  communicator = communicator_in;
104  set_owned_indices(vector_space_vector_index_set);
105  set_ghost_indices(read_write_vector_index_set);
106  }
107 
108 
109 
110  void
111  Partitioner::set_owned_indices(const IndexSet &locally_owned_indices)
112  {
113  if (Utilities::MPI::job_supports_mpi() == true)
114  {
117  }
118  else
119  {
120  my_pid = 0;
121  n_procs = 1;
122  }
123 
124  // set the local range
125  Assert(locally_owned_indices.is_contiguous() == true,
126  ExcMessage("The index set specified in locally_owned_indices "
127  "is not contiguous."));
128  locally_owned_indices.compress();
129  if (locally_owned_indices.n_elements() > 0)
131  std::pair<types::global_dof_index, types::global_dof_index>(
132  locally_owned_indices.nth_index_in_set(0),
133  locally_owned_indices.nth_index_in_set(0) +
134  locally_owned_indices.n_elements());
135  AssertThrow(
136  local_range_data.second - local_range_data.first <
137  static_cast<types::global_dof_index>(
139  ExcMessage(
140  "Index overflow: This class supports at most 2^32-1 locally owned vector entries"));
141  locally_owned_range_data.set_size(locally_owned_indices.size());
143  local_range_data.second);
145 
146  ghost_indices_data.set_size(locally_owned_indices.size());
147  }
148 
149 
150 
151  void
152  Partitioner::set_ghost_indices(const IndexSet &ghost_indices_in,
153  const IndexSet &larger_ghost_index_set)
154  {
155  // Set ghost indices from input. To be sure that no entries from the
156  // locally owned range are present, subtract the locally owned indices
157  // in any case.
158  Assert(ghost_indices_in.n_elements() == 0 ||
159  ghost_indices_in.size() == locally_owned_range_data.size(),
160  ExcDimensionMismatch(ghost_indices_in.size(),
162 
163  ghost_indices_data = ghost_indices_in;
168  AssertThrow(
170  static_cast<types::global_dof_index>(
172  ExcMessage(
173  "Index overflow: This class supports at most 2^32-1 ghost elements"));
175 
178 
179  // In the rest of this function, we determine the point-to-point
180  // communication pattern of the partitioner. We make up a list with both
181  // the processors the ghost indices actually belong to, and the indices
182  // that are locally held but ghost indices of other processors. This
183  // allows then to import and export data very easily.
184 
185  // find out the end index for each processor and communicate it (this
186  // implies the start index for the next processor)
187 #ifdef DEAL_II_WITH_MPI
188  if (n_procs < 2)
189  {
193  return;
194  }
195 
197 
198  // Allow non-zero start index for the vector. Part 1:
199  // Assume for now that the index set of rank 0 starts with 0
200  // and therefore has an increased size.
201  if (my_pid == 0)
202  my_size += local_range_data.first;
203 
204  types::global_dof_index my_shift = 0;
205  {
206  const int ierr = MPI_Exscan(&my_size,
207  &my_shift,
208  1,
210  MPI_SUM,
211  communicator);
212  AssertThrowMPI(ierr);
213  }
214 
215  // Allow non-zero start index for the vector. Part 2:
216  // We correct the assumption made above and let the
217  // index set of rank 0 actually start from the
218  // correct value, i.e. we correct the shift to
219  // its start.
220  if (my_pid == 0)
221  my_shift = local_range_data.first;
222 
223  // Fix the index start in case the index set could not give us that
224  // information.
225  if (local_range_data.first == 0 && my_shift != 0)
226  {
227  const types::global_dof_index old_local_size = local_size();
228  local_range_data.first = my_shift;
229  local_range_data.second = my_shift + old_local_size;
230  }
231 
232  std::vector<unsigned int> owning_ranks_of_ghosts(
234 
235  // set up dictionary
239  communicator,
240  owning_ranks_of_ghosts,
241  /* track origins of ghosts*/ true);
242 
243  // read dictionary by communicating with the process who owns the index
244  // in the static partition (i.e. in the dictionary). This process
245  // returns the actual owner of the index.
247  std::pair<types::global_dof_index, types::global_dof_index>,
248  unsigned int>
249  consensus_algorithm(process, communicator);
250  consensus_algorithm.run();
251 
252  {
253  ghost_targets_data = {};
254 
255  if (owning_ranks_of_ghosts.size() > 0)
256  {
257  ghost_targets_data.emplace_back(owning_ranks_of_ghosts[0], 0);
258  for (auto i : owning_ranks_of_ghosts)
259  {
260  Assert(i >= ghost_targets_data.back().first,
262  "Expect result of ConsensusAlgorithmsProcess to be "
263  "sorted"));
264  if (i == ghost_targets_data.back().first)
265  ghost_targets_data.back().second++;
266  else
267  ghost_targets_data.emplace_back(i, 1);
268  }
269  }
270  }
271 
272  // find how much the individual processes that want import from me
273  std::map<unsigned int, IndexSet> import_data = process.get_requesters();
274 
275  // count import requests and setup the compressed indices
277  import_targets_data = {};
278  import_targets_data.reserve(import_data.size());
280  import_indices_chunks_by_rank_data.reserve(import_data.size());
282  for (const auto &i : import_data)
283  if (i.second.n_elements() > 0)
284  {
285  import_targets_data.emplace_back(i.first, i.second.n_elements());
286  n_import_indices_data += i.second.n_elements();
289  i.second.n_intervals());
290  }
291 
292  // transform import indices to local index space
293  import_indices_data = {};
295  for (const auto &i : import_data)
296  {
297  Assert((i.second & locally_owned_range_data) == i.second,
298  ExcInternalError("Requested indices must be in local range"));
299  for (auto interval = i.second.begin_intervals();
300  interval != i.second.end_intervals();
301  ++interval)
302  import_indices_data.emplace_back(*interval->begin() -
303  local_range_data.first,
304  interval->last() + 1 -
305  local_range_data.first);
306  }
307 
308 # ifdef DEBUG
309 
310  // simple check: the number of processors to which we want to send
311  // ghosts and the processors to which ghosts reference should be the
312  // same
316 
317  // simple check: the number of indices to exchange should match from the
318  // ghost indices side and the import indices side
321 
322  // expensive check that the communication channel is sane -> do a ghost
323  // exchange step and see whether the ghost indices sent to us by other
324  // processes (ghost_indices) are the same as we hold locally
325  // (ghost_indices_ref).
326  std::vector<types::global_dof_index> ghost_indices_ref;
327  ghost_indices_data.fill_index_vector(ghost_indices_ref);
328  AssertDimension(ghost_indices_ref.size(), n_ghost_indices());
329  std::vector<types::global_dof_index> indices_to_send(n_import_indices());
330  std::vector<types::global_dof_index> ghost_indices(n_ghost_indices());
331  std::vector<types::global_dof_index> my_indices;
333  std::vector<MPI_Request> requests;
337  my_indices.data(), my_indices.size()),
338  make_array_view(indices_to_send),
340  requests);
342  int flag = 0;
343  const int ierr = MPI_Testall(requests.size(),
344  requests.data(),
345  &flag,
346  MPI_STATUSES_IGNORE);
347  AssertThrowMPI(ierr);
348  Assert(flag == 1,
349  ExcMessage(
350  "MPI found unfinished requests. Check communication setup"));
351 
352  for (unsigned int i = 0; i < ghost_indices.size(); ++i)
353  AssertDimension(ghost_indices[i], ghost_indices_ref[i]);
354 
355 # endif
356 
357 #endif // #ifdef DEAL_II_WITH_MPI
358 
359  if (larger_ghost_index_set.size() == 0)
360  {
362  ghost_indices_subset_data.emplace_back(local_size(),
363  local_size() +
364  n_ghost_indices());
366  }
367  else
368  {
369  AssertDimension(larger_ghost_index_set.size(),
371  Assert(
372  (larger_ghost_index_set & locally_owned_range_data).n_elements() ==
373  0,
374  ExcMessage("Ghost index set should not overlap with owned set."));
375  Assert((larger_ghost_index_set & ghost_indices_data) ==
377  ExcMessage("Larger ghost index set must contain the tight "
378  "ghost index set."));
379 
380  n_ghost_indices_in_larger_set = larger_ghost_index_set.n_elements();
381 
382  // first translate tight ghost indices into indices within the large
383  // set:
384  std::vector<unsigned int> expanded_numbering;
386  {
387  Assert(larger_ghost_index_set.is_element(index),
388  ExcMessage("The given larger ghost index set must contain "
389  "all indices in the actual index set."));
390  Assert(
391  larger_ghost_index_set.index_within_set(index) <
392  static_cast<types::global_dof_index>(
394  ExcMessage(
395  "Index overflow: This class supports at most 2^32-1 ghost elements"));
396  expanded_numbering.push_back(
397  larger_ghost_index_set.index_within_set(index));
398  }
399 
400  // now rework expanded_numbering into ranges and store in:
401  std::vector<std::pair<unsigned int, unsigned int>>
402  ghost_indices_subset;
404  ghost_targets_data.size() + 1);
405  // also populate ghost_indices_subset_chunks_by_rank_data
407  unsigned int shift = 0;
408  for (unsigned int p = 0; p < ghost_targets_data.size(); ++p)
409  {
410  unsigned int last_index = numbers::invalid_unsigned_int - 1;
411  for (unsigned int ii = 0; ii < ghost_targets_data[p].second; ii++)
412  {
413  const unsigned int i = shift + ii;
414  if (expanded_numbering[i] == last_index + 1)
415  // if contiguous, increment the end of last range:
416  ghost_indices_subset.back().second++;
417  else
418  // otherwise start a new range
419  ghost_indices_subset.emplace_back(expanded_numbering[i],
420  expanded_numbering[i] +
421  1);
422  last_index = expanded_numbering[i];
423  }
424  shift += ghost_targets_data[p].second;
426  ghost_indices_subset.size();
427  }
428  ghost_indices_subset_data = ghost_indices_subset;
429  }
430  }
431 
432 
433 
434  bool
436  {
437  // if the partitioner points to the same memory location as the calling
438  // processor
439  if (&part == this)
440  return true;
441 #ifdef DEAL_II_WITH_MPI
443  {
444  int communicators_same = 0;
445  const int ierr = MPI_Comm_compare(part.communicator,
446  communicator,
447  &communicators_same);
448  AssertThrowMPI(ierr);
449  if (!(communicators_same == MPI_IDENT ||
450  communicators_same == MPI_CONGRUENT))
451  return false;
452  }
453 #endif
454  return (global_size == part.global_size &&
457  }
458 
459 
460 
461  bool
463  {
464  return Utilities::MPI::min(static_cast<int>(is_compatible(part)),
465  communicator) == 1;
466  }
467 
468 
469 
470  std::size_t
472  {
473  std::size_t memory = (3 * sizeof(types::global_dof_index) +
474  4 * sizeof(unsigned int) + sizeof(MPI_Comm));
483  memory +=
486  return memory;
487  }
488 
489  } // namespace MPI
490 
491 } // end of namespace Utilities
492 
493 
494 
495 // explicit instantiations from .templates.h file
496 #include "partitioner.inst"
497 
Utilities::MPI::Partitioner::n_ghost_indices_in_larger_set
unsigned int n_ghost_indices_in_larger_set
Definition: partitioner.h:728
Utilities::MPI::Partitioner::set_ghost_indices
void set_ghost_indices(const IndexSet &ghost_indices, const IndexSet &larger_ghost_index_set=IndexSet())
Definition: partitioner.cc:152
IndexSet::compress
void compress() const
Definition: index_set.h:1643
IndexSet
Definition: index_set.h:74
IndexSet::subtract_set
void subtract_set(const IndexSet &other)
Definition: index_set.cc:258
Utilities::MPI::Partitioner::Partitioner
Partitioner()
Definition: partitioner.cc:26
Utilities::MPI::Partitioner
Definition: partitioner.h:195
DEAL_II_DOF_INDEX_MPI_TYPE
#define DEAL_II_DOF_INDEX_MPI_TYPE
Definition: types.h:86
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::get_requesters
std::map< unsigned int, IndexSet > get_requesters()
Definition: mpi_compute_index_owner_internal.h:782
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
Definition: mpi_compute_index_owner_internal.h:538
IndexSet::is_contiguous
bool is_contiguous() const
Definition: index_set.h:1816
Utilities::MPI::Partitioner::import_targets_data
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
Definition: partitioner.h:716
ArrayView
Definition: array_view.h:77
IndexSet::size
size_type size() const
Definition: index_set.h:1635
partitioner.h
MPI_Comm
make_array_view
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition: array_view.h:607
Utilities::MPI::Partitioner::global_size
types::global_dof_index global_size
Definition: partitioner.h:652
Utilities::MPI::Partitioner::size
types::global_dof_index size() const
Utilities::MPI::Partitioner::n_import_indices_data
unsigned int n_import_indices_data
Definition: partitioner.h:710
Utilities::MPI::ConsensusAlgorithms::Selector
Definition: mpi_consensus_algorithms.h:464
IndexSet::n_elements
size_type n_elements() const
Definition: index_set.h:1833
IndexSet::set_size
void set_size(const size_type size)
Definition: index_set.h:1623
Utilities::MPI::Partitioner::import_indices_data
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
Definition: partitioner.h:690
Utilities::MPI::Partitioner::ghost_targets_data
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
Definition: partitioner.h:682
Utilities::MPI::Partitioner::have_ghost_indices
bool have_ghost_indices
Definition: partitioner.h:762
Utilities::MPI::Partitioner::export_to_ghosted_array_finish
void export_to_ghosted_array_finish(const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
Utilities::MPI::Partitioner::n_ghost_indices_data
unsigned int n_ghost_indices_data
Definition: partitioner.h:676
Utilities::MPI::Partitioner::set_owned_indices
void set_owned_indices(const IndexSet &locally_owned_indices)
Definition: partitioner.cc:111
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
Utilities::MPI::Partitioner::ghost_indices_subset_data
std::vector< std::pair< unsigned int, unsigned int > > ghost_indices_subset_data
Definition: partitioner.h:742
Utilities::MPI::Partitioner::is_compatible
bool is_compatible(const Partitioner &part) const
Definition: partitioner.cc:435
Utilities::MPI::Partitioner::memory_consumption
std::size_t memory_consumption() const
Definition: partitioner.cc:471
Utilities::MPI::Partitioner::n_ghost_indices
unsigned int n_ghost_indices() const
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
Utilities::MPI::Partitioner::ghost_indices_subset_chunks_by_rank_data
std::vector< unsigned int > ghost_indices_subset_chunks_by_rank_data
Definition: partitioner.h:734
IndexSet::is_element
bool is_element(const size_type index) const
Definition: index_set.h:1766
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
Utilities::MPI::Partitioner::locally_owned_range_data
IndexSet locally_owned_range_data
Definition: partitioner.h:657
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
Utilities::MPI::Partitioner::n_import_indices
unsigned int n_import_indices() const
types
Definition: types.h:31
Utilities::MPI::Partitioner::local_size
unsigned int local_size() const
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Utilities::MPI::Partitioner::communicator
MPI_Comm communicator
Definition: partitioner.h:757
Utilities::MPI::Partitioner::is_globally_compatible
bool is_globally_compatible(const Partitioner &part) const
Definition: partitioner.cc:462
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Utilities::MPI::this_mpi_process
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
IndexSet::nth_index_in_set
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1881
Utilities::MPI::Partitioner::ghost_indices_data
IndexSet ghost_indices_data
Definition: partitioner.h:670
GridTools::shift
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:817
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
Utilities::MPI::Partitioner::export_to_ghosted_array_start
void export_to_ghosted_array_start(const unsigned int communication_channel, const ArrayView< const Number, MemorySpaceType > &locally_owned_array, const ArrayView< Number, MemorySpaceType > &temporary_storage, const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Utilities::MPI::Partitioner::reinit
virtual void reinit(const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator) override
Definition: partitioner.cc:98
IndexSet::index_within_set
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1922
Utilities::MPI::Partitioner::n_procs
unsigned int n_procs
Definition: partitioner.h:752
Utilities::MPI::Partitioner::ghost_indices
const IndexSet & ghost_indices() const
Utilities::MPI::Partitioner::my_pid
unsigned int my_pid
Definition: partitioner.h:747
IndexSet::fill_index_vector
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:507
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Utilities::MPI::job_supports_mpi
bool job_supports_mpi()
Definition: mpi.cc:1030
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
Utilities
Definition: cuda.h:31
IndexSet::add_range
void add_range(const size_type begin, const size_type end)
Definition: index_set.h:1674
Utilities::MPI::Partitioner::import_indices_chunks_by_rank_data
std::vector< unsigned int > import_indices_chunks_by_rank_data
Definition: partitioner.h:722
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
mpi_compute_index_owner_internal.h
Utilities::MPI::Partitioner::local_range_data
std::pair< types::global_dof_index, types::global_dof_index > local_range_data
Definition: partitioner.h:664