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