Reference documentation for deal.II version 9.3.3
\(\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 - 2021 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
22namespace Utilities
23{
24 namespace MPI
25 {
27 : global_size(0)
28 , local_range_data(
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(
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 {
54 locally_owned_range_data.add_range(0, size);
55 locally_owned_range_data.compress();
56 ghost_indices_data.set_size(size);
57 }
58
59
60
62 const types::global_dof_index ghost_size,
63 const MPI_Comm & communicator)
64 : global_size(Utilities::MPI::sum<types::global_dof_index>(local_size,
65 communicator))
66 , locally_owned_range_data(global_size)
67 , local_range_data{0, local_size}
68 , n_ghost_indices_data(ghost_size)
69 , n_import_indices_data(0)
70 , n_ghost_indices_in_larger_set(0)
71 , my_pid(Utilities::MPI::this_mpi_process(communicator))
72 , n_procs(Utilities::MPI::n_mpi_processes(communicator))
73 , communicator(communicator)
74 , have_ghost_indices(true)
75 {
76 types::global_dof_index prefix_sum = 0;
77
78#ifdef DEAL_II_WITH_MPI
79 const int ierr =
80 MPI_Exscan(&local_size,
81 &prefix_sum,
82 1,
83 Utilities::MPI::internal::mpi_type_id(&prefix_sum),
84 MPI_SUM,
85 communicator);
86 AssertThrowMPI(ierr);
87#endif
88
89 local_range_data = {prefix_sum, prefix_sum + local_size};
90
91 locally_owned_range_data.add_range(prefix_sum, prefix_sum + local_size);
92 locally_owned_range_data.compress();
93 }
94
95
96
97 Partitioner::Partitioner(const IndexSet &locally_owned_indices,
98 const IndexSet &ghost_indices_in,
99 const MPI_Comm &communicator_in)
100 : global_size(
101 static_cast<types::global_dof_index>(locally_owned_indices.size()))
102 , n_ghost_indices_data(0)
103 , n_import_indices_data(0)
104 , n_ghost_indices_in_larger_set(0)
105 , my_pid(0)
106 , n_procs(1)
107 , communicator(communicator_in)
108 , have_ghost_indices(false)
109 {
110 set_owned_indices(locally_owned_indices);
111 set_ghost_indices(ghost_indices_in);
112 }
113
114
115
116 Partitioner::Partitioner(const IndexSet &locally_owned_indices,
117 const MPI_Comm &communicator_in)
118 : global_size(
119 static_cast<types::global_dof_index>(locally_owned_indices.size()))
120 , n_ghost_indices_data(0)
121 , n_import_indices_data(0)
122 , n_ghost_indices_in_larger_set(0)
123 , my_pid(0)
124 , n_procs(1)
125 , communicator(communicator_in)
126 , have_ghost_indices(false)
127 {
128 set_owned_indices(locally_owned_indices);
129 }
130
131
132
133 void
134 Partitioner::reinit(const IndexSet &vector_space_vector_index_set,
135 const IndexSet &read_write_vector_index_set,
136 const MPI_Comm &communicator_in)
137 {
138 have_ghost_indices = false;
139 communicator = communicator_in;
140 set_owned_indices(vector_space_vector_index_set);
141 set_ghost_indices(read_write_vector_index_set);
142 }
143
144
145
146 void
147 Partitioner::set_owned_indices(const IndexSet &locally_owned_indices)
148 {
150 {
151 my_pid = Utilities::MPI::this_mpi_process(communicator);
152 n_procs = Utilities::MPI::n_mpi_processes(communicator);
153 }
154 else
155 {
156 my_pid = 0;
157 n_procs = 1;
158 }
159
160 // set the local range
161 Assert(locally_owned_indices.is_contiguous() == true,
162 ExcMessage("The index set specified in locally_owned_indices "
163 "is not contiguous."));
164 locally_owned_indices.compress();
165 if (locally_owned_indices.n_elements() > 0)
166 local_range_data =
167 std::pair<types::global_dof_index, types::global_dof_index>(
168 locally_owned_indices.nth_index_in_set(0),
169 locally_owned_indices.nth_index_in_set(0) +
170 locally_owned_indices.n_elements());
172 local_range_data.second - local_range_data.first <
173 static_cast<types::global_dof_index>(
176 "Index overflow: This class supports at most 2^32-1 locally owned vector entries"));
177 locally_owned_range_data.set_size(locally_owned_indices.size());
178 locally_owned_range_data.add_range(local_range_data.first,
179 local_range_data.second);
180 locally_owned_range_data.compress();
181
182 ghost_indices_data.set_size(locally_owned_indices.size());
183 }
184
185
186
187 void
188 Partitioner::set_ghost_indices(const IndexSet &ghost_indices_in,
189 const IndexSet &larger_ghost_index_set)
190 {
191 // Set ghost indices from input. To be sure that no entries from the
192 // locally owned range are present, subtract the locally owned indices
193 // in any case.
194 Assert(ghost_indices_in.n_elements() == 0 ||
195 ghost_indices_in.size() == locally_owned_range_data.size(),
196 ExcDimensionMismatch(ghost_indices_in.size(),
197 locally_owned_range_data.size()));
198
199 ghost_indices_data = ghost_indices_in;
200 if (ghost_indices_data.size() != locally_owned_range_data.size())
201 ghost_indices_data.set_size(locally_owned_range_data.size());
202 ghost_indices_data.subtract_set(locally_owned_range_data);
203 ghost_indices_data.compress();
205 ghost_indices_data.n_elements() <
206 static_cast<types::global_dof_index>(
209 "Index overflow: This class supports at most 2^32-1 ghost elements"));
210 n_ghost_indices_data = ghost_indices_data.n_elements();
211
212 have_ghost_indices =
213 Utilities::MPI::max(n_ghost_indices_data, communicator) > 0;
214
215 // In the rest of this function, we determine the point-to-point
216 // communication pattern of the partitioner. We make up a list with both
217 // the processors the ghost indices actually belong to, and the indices
218 // that are locally held but ghost indices of other processors. This
219 // allows then to import and export data very easily.
220
221 // find out the end index for each processor and communicate it (this
222 // implies the start index for the next processor)
223#ifdef DEAL_II_WITH_MPI
224 if (n_procs < 2)
225 {
226 Assert(ghost_indices_data.n_elements() == 0, ExcInternalError());
227 Assert(n_import_indices_data == 0, ExcInternalError());
228 Assert(n_ghost_indices_data == 0, ExcInternalError());
229 return;
230 }
231
232 types::global_dof_index my_size = locally_owned_size();
233
234 // Allow non-zero start index for the vector. Part 1:
235 // Assume for now that the index set of rank 0 starts with 0
236 // and therefore has an increased size.
237 if (my_pid == 0)
238 my_size += local_range_data.first;
239
240 types::global_dof_index my_shift = 0;
241 {
242 const int ierr = MPI_Exscan(&my_size,
243 &my_shift,
244 1,
246 MPI_SUM,
247 communicator);
248 AssertThrowMPI(ierr);
249 }
250
251 // Allow non-zero start index for the vector. Part 2:
252 // We correct the assumption made above and let the
253 // index set of rank 0 actually start from the
254 // correct value, i.e. we correct the shift to
255 // its start.
256 if (my_pid == 0)
257 my_shift = local_range_data.first;
258
259 // Fix the index start in case the index set could not give us that
260 // information.
261 if (local_range_data.first == 0 && my_shift != 0)
262 {
263 const types::global_dof_index old_locally_owned_size =
264 locally_owned_size();
265 local_range_data.first = my_shift;
266 local_range_data.second = my_shift + old_locally_owned_size;
267 }
268
269 std::vector<unsigned int> owning_ranks_of_ghosts(
270 ghost_indices_data.n_elements());
271
272 // set up dictionary
273 internal::ComputeIndexOwner::ConsensusAlgorithmsPayload process(
274 locally_owned_range_data,
275 ghost_indices_data,
276 communicator,
277 owning_ranks_of_ghosts,
278 /* track origins of ghosts*/ true);
279
280 // read dictionary by communicating with the process who owns the index
281 // in the static partition (i.e. in the dictionary). This process
282 // returns the actual owner of the index.
283 ConsensusAlgorithms::Selector<
284 std::pair<types::global_dof_index, types::global_dof_index>,
285 unsigned int>
286 consensus_algorithm(process, communicator);
287 consensus_algorithm.run();
288
289 {
290 ghost_targets_data = {};
291
292 if (owning_ranks_of_ghosts.size() > 0)
293 {
294 ghost_targets_data.emplace_back(owning_ranks_of_ghosts[0], 0);
295 for (auto i : owning_ranks_of_ghosts)
296 {
297 Assert(i >= ghost_targets_data.back().first,
299 "Expect result of ConsensusAlgorithmsProcess to be "
300 "sorted"));
301 if (i == ghost_targets_data.back().first)
302 ghost_targets_data.back().second++;
303 else
304 ghost_targets_data.emplace_back(i, 1);
305 }
306 }
307 }
308
309 // find how much the individual processes that want import from me
310 std::map<unsigned int, IndexSet> import_data = process.get_requesters();
311
312 // count import requests and setup the compressed indices
313 n_import_indices_data = 0;
314 import_targets_data = {};
315 import_targets_data.reserve(import_data.size());
316 import_indices_chunks_by_rank_data = {};
317 import_indices_chunks_by_rank_data.reserve(import_data.size());
318 import_indices_chunks_by_rank_data.resize(1);
319 for (const auto &i : import_data)
320 if (i.second.n_elements() > 0)
321 {
322 import_targets_data.emplace_back(i.first, i.second.n_elements());
323 n_import_indices_data += i.second.n_elements();
324 import_indices_chunks_by_rank_data.push_back(
325 import_indices_chunks_by_rank_data.back() +
326 i.second.n_intervals());
327 }
328
329 // transform import indices to local index space
330 import_indices_data = {};
331 import_indices_data.reserve(import_indices_chunks_by_rank_data.back());
332 for (const auto &i : import_data)
333 {
334 Assert((i.second & locally_owned_range_data) == i.second,
335 ExcInternalError("Requested indices must be in local range"));
336 for (auto interval = i.second.begin_intervals();
337 interval != i.second.end_intervals();
338 ++interval)
339 import_indices_data.emplace_back(*interval->begin() -
340 local_range_data.first,
341 interval->last() + 1 -
342 local_range_data.first);
343 }
344
345# ifdef DEBUG
346
347 // simple check: the number of processors to which we want to send
348 // ghosts and the processors to which ghosts reference should be the
349 // same
351 Utilities::MPI::sum(import_targets_data.size(), communicator),
352 Utilities::MPI::sum(ghost_targets_data.size(), communicator));
353
354 // simple check: the number of indices to exchange should match from the
355 // ghost indices side and the import indices side
356 AssertDimension(Utilities::MPI::sum(n_import_indices_data, communicator),
357 Utilities::MPI::sum(n_ghost_indices_data, communicator));
358
359 // expensive check that the communication channel is sane -> do a ghost
360 // exchange step and see whether the ghost indices sent to us by other
361 // processes (ghost_indices) are the same as we hold locally
362 // (ghost_indices_ref).
363 std::vector<types::global_dof_index> ghost_indices_ref;
364 ghost_indices_data.fill_index_vector(ghost_indices_ref);
365 AssertDimension(ghost_indices_ref.size(), n_ghost_indices());
366 std::vector<types::global_dof_index> indices_to_send(n_import_indices());
367 std::vector<types::global_dof_index> ghost_indices(n_ghost_indices());
368 std::vector<types::global_dof_index> my_indices;
369 locally_owned_range_data.fill_index_vector(my_indices);
370 std::vector<MPI_Request> requests;
371 n_ghost_indices_in_larger_set = n_ghost_indices_data;
372 export_to_ghosted_array_start(127,
374 my_indices.data(), my_indices.size()),
375 make_array_view(indices_to_send),
376 make_array_view(ghost_indices),
377 requests);
378 export_to_ghosted_array_finish(make_array_view(ghost_indices), requests);
379 int flag = 0;
380 const int ierr = MPI_Testall(requests.size(),
381 requests.data(),
382 &flag,
383 MPI_STATUSES_IGNORE);
384 AssertThrowMPI(ierr);
385 Assert(flag == 1,
387 "MPI found unfinished requests. Check communication setup"));
388
389 for (unsigned int i = 0; i < ghost_indices.size(); ++i)
390 AssertDimension(ghost_indices[i], ghost_indices_ref[i]);
391
392# endif
393
394#endif // #ifdef DEAL_II_WITH_MPI
395
396 if (larger_ghost_index_set.size() == 0)
397 {
398 ghost_indices_subset_chunks_by_rank_data.clear();
399 ghost_indices_subset_data.emplace_back(0, n_ghost_indices());
400 n_ghost_indices_in_larger_set = n_ghost_indices_data;
401 }
402 else
403 {
404 AssertDimension(larger_ghost_index_set.size(),
405 ghost_indices_data.size());
406 Assert(
407 (larger_ghost_index_set & locally_owned_range_data).n_elements() ==
408 0,
409 ExcMessage("Ghost index set should not overlap with owned set."));
410 Assert((larger_ghost_index_set & ghost_indices_data) ==
411 ghost_indices_data,
412 ExcMessage("Larger ghost index set must contain the tight "
413 "ghost index set."));
414
415 n_ghost_indices_in_larger_set = larger_ghost_index_set.n_elements();
416
417 // first translate tight ghost indices into indices within the large
418 // set:
419 std::vector<unsigned int> expanded_numbering;
420 for (::IndexSet::size_type index : ghost_indices_data)
421 {
422 Assert(larger_ghost_index_set.is_element(index),
423 ExcMessage("The given larger ghost index set must contain "
424 "all indices in the actual index set."));
425 Assert(
426 larger_ghost_index_set.index_within_set(index) <
427 static_cast<types::global_dof_index>(
430 "Index overflow: This class supports at most 2^32-1 ghost elements"));
431 expanded_numbering.push_back(
432 larger_ghost_index_set.index_within_set(index));
433 }
434
435 // now rework expanded_numbering into ranges and store in:
436 std::vector<std::pair<unsigned int, unsigned int>>
437 ghost_indices_subset;
438 ghost_indices_subset_chunks_by_rank_data.resize(
439 ghost_targets_data.size() + 1);
440 // also populate ghost_indices_subset_chunks_by_rank_data
441 ghost_indices_subset_chunks_by_rank_data[0] = 0;
442 unsigned int shift = 0;
443 for (unsigned int p = 0; p < ghost_targets_data.size(); ++p)
444 {
445 unsigned int last_index = numbers::invalid_unsigned_int - 1;
446 for (unsigned int ii = 0; ii < ghost_targets_data[p].second; ii++)
447 {
448 const unsigned int i = shift + ii;
449 if (expanded_numbering[i] == last_index + 1)
450 // if contiguous, increment the end of last range:
451 ghost_indices_subset.back().second++;
452 else
453 // otherwise start a new range
454 ghost_indices_subset.emplace_back(expanded_numbering[i],
455 expanded_numbering[i] +
456 1);
457 last_index = expanded_numbering[i];
458 }
459 shift += ghost_targets_data[p].second;
460 ghost_indices_subset_chunks_by_rank_data[p + 1] =
461 ghost_indices_subset.size();
462 }
463 ghost_indices_subset_data = ghost_indices_subset;
464 }
465 }
466
467
468
469 bool
470 Partitioner::is_compatible(const Partitioner &part) const
471 {
472 // if the partitioner points to the same memory location as the calling
473 // processor
474 if (&part == this)
475 return true;
476#ifdef DEAL_II_WITH_MPI
478 {
479 int communicators_same = 0;
480 const int ierr = MPI_Comm_compare(part.communicator,
481 communicator,
482 &communicators_same);
483 AssertThrowMPI(ierr);
484 if (!(communicators_same == MPI_IDENT ||
485 communicators_same == MPI_CONGRUENT))
486 return false;
487 }
488#endif
489 return (global_size == part.global_size &&
490 local_range_data == part.local_range_data &&
491 ghost_indices_data == part.ghost_indices_data);
492 }
493
494
495
496 bool
497 Partitioner::is_globally_compatible(const Partitioner &part) const
498 {
499 return Utilities::MPI::min(static_cast<int>(is_compatible(part)),
500 communicator) == 1;
501 }
502
503
504
505 std::size_t
507 {
508 std::size_t memory = (3 * sizeof(types::global_dof_index) +
509 4 * sizeof(unsigned int) + sizeof(MPI_Comm));
510 memory += MemoryConsumption::memory_consumption(locally_owned_range_data);
511 memory += MemoryConsumption::memory_consumption(ghost_targets_data);
512 memory += MemoryConsumption::memory_consumption(import_targets_data);
513 memory += MemoryConsumption::memory_consumption(import_indices_data);
515 import_indices_chunks_by_rank_data);
517 ghost_indices_subset_chunks_by_rank_data);
518 memory +=
519 MemoryConsumption::memory_consumption(ghost_indices_subset_data);
520 memory += MemoryConsumption::memory_consumption(ghost_indices_data);
521 return memory;
522 }
523
524 } // namespace MPI
525
526} // end of namespace Utilities
527
528
529
530// explicit instantiations from .templates.h file
531#include "partitioner.inst"
532
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:697
bool is_contiguous() const
Definition: index_set.h:1815
size_type size() const
Definition: index_set.h:1634
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1921
size_type n_elements() const
Definition: index_set.h:1832
bool is_element(const size_type index) const
Definition: index_set.h:1765
void set_size(const size_type size)
Definition: index_set.h:1622
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1880
void compress() const
Definition: index_set.h:1642
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1746
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:2022
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
T min(const T &t, const MPI_Comm &mpi_communicator)
bool job_supports_mpi()
Definition: mpi.cc:1097
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
T max(const T &t, const MPI_Comm &mpi_communicator)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
static const unsigned int invalid_unsigned_int
Definition: types.h:196
STL namespace.
Definition: types.h:32
unsigned int global_dof_index
Definition: types.h:76
#define DEAL_II_DOF_INDEX_MPI_TYPE
Definition: types.h:86