18#include <deal.II/base/partitioner.templates.h>
30 , n_ghost_indices_data(0)
31 , n_import_indices_data(0)
32 , n_ghost_indices_in_larger_set(0)
35 , communicator(MPI_COMM_SELF)
36 , have_ghost_indices(false)
43 , locally_owned_range_data(size)
46 , n_ghost_indices_data(0)
47 , n_import_indices_data(0)
48 , n_ghost_indices_in_larger_set(0)
51 , communicator(MPI_COMM_SELF)
52 , have_ghost_indices(false)
54 locally_owned_range_data.add_range(0, size);
55 locally_owned_range_data.compress();
56 ghost_indices_data.set_size(size);
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)
73 , communicator(communicator)
74 , have_ghost_indices(true)
78#ifdef DEAL_II_WITH_MPI
80 MPI_Exscan(&local_size,
83 Utilities::MPI::internal::mpi_type_id(&prefix_sum),
89 local_range_data = {prefix_sum, prefix_sum + local_size};
91 locally_owned_range_data.add_range(prefix_sum, prefix_sum + local_size);
92 locally_owned_range_data.compress();
102 , n_ghost_indices_data(0)
103 , n_import_indices_data(0)
104 , n_ghost_indices_in_larger_set(0)
107 , communicator(communicator_in)
108 , have_ghost_indices(false)
110 set_owned_indices(locally_owned_indices);
111 set_ghost_indices(ghost_indices_in);
120 , n_ghost_indices_data(0)
121 , n_import_indices_data(0)
122 , n_ghost_indices_in_larger_set(0)
125 , communicator(communicator_in)
126 , have_ghost_indices(false)
128 set_owned_indices(locally_owned_indices);
135 const IndexSet &read_write_vector_index_set,
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);
147 Partitioner::set_owned_indices(
const IndexSet &locally_owned_indices)
162 ExcMessage(
"The index set specified in locally_owned_indices "
163 "is not contiguous."));
167 std::pair<types::global_dof_index, types::global_dof_index>(
172 local_range_data.second - local_range_data.first <
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();
182 ghost_indices_data.set_size(locally_owned_indices.
size());
188 Partitioner::set_ghost_indices(
const IndexSet &ghost_indices_in,
189 const IndexSet &larger_ghost_index_set)
195 ghost_indices_in.
size() == locally_owned_range_data.size(),
197 locally_owned_range_data.size()));
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() <
209 "Index overflow: This class supports at most 2^32-1 ghost elements"));
210 n_ghost_indices_data = ghost_indices_data.n_elements();
223#ifdef DEAL_II_WITH_MPI
238 my_size += local_range_data.first;
242 const int ierr = MPI_Exscan(&my_size,
257 my_shift = local_range_data.first;
261 if (local_range_data.first == 0 && my_shift != 0)
264 locally_owned_size();
265 local_range_data.first = my_shift;
266 local_range_data.second = my_shift + old_locally_owned_size;
269 std::vector<unsigned int> owning_ranks_of_ghosts(
270 ghost_indices_data.n_elements());
273 internal::ComputeIndexOwner::ConsensusAlgorithmsPayload process(
274 locally_owned_range_data,
277 owning_ranks_of_ghosts,
283 ConsensusAlgorithms::Selector<
284 std::pair<types::global_dof_index, types::global_dof_index>,
286 consensus_algorithm(process, communicator);
287 consensus_algorithm.run();
290 ghost_targets_data = {};
292 if (owning_ranks_of_ghosts.size() > 0)
294 ghost_targets_data.emplace_back(owning_ranks_of_ghosts[0], 0);
295 for (
auto i : owning_ranks_of_ghosts)
297 Assert(i >= ghost_targets_data.back().first,
299 "Expect result of ConsensusAlgorithmsProcess to be "
301 if (i == ghost_targets_data.back().first)
302 ghost_targets_data.back().second++;
304 ghost_targets_data.emplace_back(i, 1);
310 std::map<unsigned int, IndexSet> import_data = process.get_requesters();
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)
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());
330 import_indices_data = {};
331 import_indices_data.reserve(import_indices_chunks_by_rank_data.back());
332 for (
const auto &i : import_data)
334 Assert((i.second & locally_owned_range_data) == i.second,
336 for (
auto interval = i.second.begin_intervals();
337 interval != i.second.end_intervals();
339 import_indices_data.emplace_back(*interval->begin() -
340 local_range_data.first,
341 interval->last() + 1 -
342 local_range_data.first);
363 std::vector<types::global_dof_index> ghost_indices_ref;
364 ghost_indices_data.fill_index_vector(ghost_indices_ref);
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()),
378 export_to_ghosted_array_finish(
make_array_view(ghost_indices), requests);
380 const int ierr = MPI_Testall(requests.size(),
383 MPI_STATUSES_IGNORE);
387 "MPI found unfinished requests. Check communication setup"));
389 for (
unsigned int i = 0; i < ghost_indices.size(); ++i)
396 if (larger_ghost_index_set.
size() == 0)
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;
405 ghost_indices_data.size());
407 (larger_ghost_index_set & locally_owned_range_data).n_elements() ==
409 ExcMessage(
"Ghost index set should not overlap with owned set."));
410 Assert((larger_ghost_index_set & ghost_indices_data) ==
412 ExcMessage(
"Larger ghost index set must contain the tight "
413 "ghost index set."));
415 n_ghost_indices_in_larger_set = larger_ghost_index_set.
n_elements();
419 std::vector<unsigned int> expanded_numbering;
423 ExcMessage(
"The given larger ghost index set must contain "
424 "all indices in the actual index set."));
430 "Index overflow: This class supports at most 2^32-1 ghost elements"));
431 expanded_numbering.push_back(
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);
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)
446 for (
unsigned int ii = 0; ii < ghost_targets_data[p].second; ii++)
448 const unsigned int i =
shift + ii;
449 if (expanded_numbering[i] == last_index + 1)
451 ghost_indices_subset.back().second++;
454 ghost_indices_subset.emplace_back(expanded_numbering[i],
455 expanded_numbering[i] +
457 last_index = expanded_numbering[i];
459 shift += ghost_targets_data[p].second;
460 ghost_indices_subset_chunks_by_rank_data[p + 1] =
461 ghost_indices_subset.size();
463 ghost_indices_subset_data = ghost_indices_subset;
470 Partitioner::is_compatible(
const Partitioner &part)
const
476#ifdef DEAL_II_WITH_MPI
479 int communicators_same = 0;
480 const int ierr = MPI_Comm_compare(part.communicator,
482 &communicators_same);
484 if (!(communicators_same == MPI_IDENT ||
485 communicators_same == MPI_CONGRUENT))
489 return (global_size == part.global_size &&
490 local_range_data == part.local_range_data &&
491 ghost_indices_data == part.ghost_indices_data);
497 Partitioner::is_globally_compatible(
const Partitioner &part)
const
509 4 *
sizeof(
unsigned int) +
sizeof(
MPI_Comm));
515 import_indices_chunks_by_rank_data);
517 ghost_indices_subset_chunks_by_rank_data);
531#include "partitioner.inst"
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
bool is_contiguous() const
size_type index_within_set(const size_type global_index) const
size_type n_elements() const
bool is_element(const size_type index) const
void set_size(const size_type size)
size_type nth_index_in_set(const size_type local_index) const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
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)
T min(const T &t, const MPI_Comm &mpi_communicator)
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
T max(const T &t, const MPI_Comm &mpi_communicator)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static const unsigned int invalid_unsigned_int
unsigned int global_dof_index
#define DEAL_II_DOF_INDEX_MPI_TYPE