18#include <deal.II/base/partitioner.templates.h>
20#include <boost/serialization/utility.hpp>
33 std::pair<
types::global_dof_index,
types::global_dof_index>(0, 0))
34 , n_ghost_indices_data(0)
35 , n_import_indices_data(0)
36 , n_ghost_indices_in_larger_set(0)
39 , communicator(MPI_COMM_SELF)
40 , have_ghost_indices(false)
47 , locally_owned_range_data(size)
49 std::pair<
types::global_dof_index,
types::global_dof_index>(0, size))
50 , n_ghost_indices_data(0)
51 , n_import_indices_data(0)
52 , n_ghost_indices_in_larger_set(0)
55 , communicator(MPI_COMM_SELF)
56 , have_ghost_indices(false)
70 , locally_owned_range_data(global_size)
71 , local_range_data{0, local_size}
72 , n_ghost_indices_data(ghost_size)
73 , n_import_indices_data(0)
74 , n_ghost_indices_in_larger_set(0)
77 , communicator(communicator)
78 , have_ghost_indices(true)
82#ifdef DEAL_II_WITH_MPI
105 static_cast<
types::global_dof_index>(locally_owned_indices.size()))
106 , n_ghost_indices_data(0)
107 , n_import_indices_data(0)
108 , n_ghost_indices_in_larger_set(0)
111 , communicator(communicator_in)
112 , have_ghost_indices(false)
123 static_cast<
types::global_dof_index>(locally_owned_indices.size()))
124 , n_ghost_indices_data(0)
125 , n_import_indices_data(0)
126 , n_ghost_indices_in_larger_set(0)
129 , communicator(communicator_in)
130 , have_ghost_indices(false)
139 const IndexSet &read_write_vector_index_set,
166 ExcMessage(
"The index set specified in locally_owned_indices "
167 "is not contiguous."));
171 std::pair<types::global_dof_index, types::global_dof_index>(
178 std::numeric_limits<unsigned int>::max()),
180 "Index overflow: This class supports at most 2^32-1 locally owned vector entries"));
193 const IndexSet &larger_ghost_index_set)
211 std::numeric_limits<unsigned int>::max()),
213 "Index overflow: This class supports at most 2^32-1 ghost elements"));
227#ifdef DEAL_II_WITH_MPI
246 const int ierr = MPI_Exscan(&my_size,
273 std::vector<unsigned int> owning_ranks_of_ghosts(
281 owning_ranks_of_ghosts,
289 std::pair<types::global_dof_index, types::global_dof_index>>,
290 std::vector<unsigned int>>
297 if (owning_ranks_of_ghosts.size() > 0)
300 for (
auto i : owning_ranks_of_ghosts)
304 "Expect result of ConsensusAlgorithms::Process to be "
315 std::map<unsigned int, IndexSet> import_data = process.
get_requesters();
324 for (
const auto &i : import_data)
325 if (i.second.n_elements() > 0)
331 i.second.n_intervals());
337 for (
const auto &i : import_data)
341 for (
auto interval = i.second.begin_intervals();
342 interval != i.second.end_intervals();
346 interval->last() + 1 -
368 const std::vector<types::global_dof_index> ghost_indices_ref =
374 const std::vector<types::global_dof_index> my_indices =
376 std::vector<MPI_Request> requests;
380 my_indices.data(), my_indices.size()),
386 const int ierr = MPI_Testall(requests.size(),
389 MPI_STATUSES_IGNORE);
393 "MPI found unfinished requests. Check communication setup"));
402 if (larger_ghost_index_set.
size() == 0)
415 ExcMessage(
"Ghost index set should not overlap with owned set."));
418 ExcMessage(
"Larger ghost index set must contain the tight "
419 "ghost index set."));
425 std::vector<unsigned int> expanded_numbering;
429 ExcMessage(
"The given larger ghost index set must contain "
430 "all indices in the actual index set."));
434 std::numeric_limits<unsigned int>::max()),
436 "Index overflow: This class supports at most 2^32-1 ghost elements"));
437 expanded_numbering.push_back(
442 std::vector<std::pair<unsigned int, unsigned int>>
443 ghost_indices_subset;
448 unsigned int shift = 0;
454 const unsigned int i = shift + ii;
455 if (expanded_numbering[i] == last_index + 1)
457 ghost_indices_subset.back().second++;
460 ghost_indices_subset.emplace_back(expanded_numbering[i],
461 expanded_numbering[i] +
463 last_index = expanded_numbering[i];
467 ghost_indices_subset.size();
482#ifdef DEAL_II_WITH_MPI
485 int communicators_same = 0;
488 &communicators_same);
490 if (!(communicators_same == MPI_IDENT ||
491 communicators_same == MPI_CONGRUENT))
548 for (
unsigned int i = 0; i < n_import_targets; ++i)
551 std::vector<std::pair<unsigned int, unsigned int>>::const_iterator
556 std::vector<unsigned int> import_indices_plain_host;
557 for (; my_imports != end_my_imports; ++my_imports)
559 const unsigned int chunk_size =
560 my_imports->second - my_imports->first;
561 for (
unsigned int j = 0; j < chunk_size; ++j)
562 import_indices_plain_host.push_back(my_imports->first + j);
566 const auto chunk_size = import_indices_plain_host.size();
571 Kokkos::View<unsigned int *, Kokkos::HostSpace>(
572 import_indices_plain_host.data(), chunk_size));
583#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)
void subtract_set(const IndexSet &other)
void add_range(const size_type begin, const size_type end)
size_type nth_index_in_set(const size_type local_index) const
std::size_t memory_consumption() const
std::vector< size_type > get_index_vector() const
virtual std::vector< unsigned int > run(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm) override
unsigned int n_ghost_indices_data
std::size_t memory_consumption() const
std::vector< Kokkos::View< unsigned int *, MemorySpace::Default::kokkos_space > > import_indices_plain_dev
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
bool is_globally_compatible(const Partitioner &part) const
void set_owned_indices(const IndexSet &locally_owned_indices)
IndexSet ghost_indices_data
unsigned int n_ghost_indices_in_larger_set
unsigned int n_ghost_indices() const
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
unsigned int local_size() const
const IndexSet & ghost_indices() const
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
virtual void reinit(const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm communicator) override
void initialize_import_indices_plain_dev() const
types::global_dof_index global_size
unsigned int locally_owned_size() const
std::pair< types::global_dof_index, types::global_dof_index > local_range_data
std::vector< unsigned int > import_indices_chunks_by_rank_data
unsigned int n_import_indices_data
void set_ghost_indices(const IndexSet &ghost_indices, const IndexSet &larger_ghost_index_set=IndexSet())
void export_to_ghosted_array_finish(const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
unsigned int n_import_indices() const
std::vector< unsigned int > ghost_indices_subset_chunks_by_rank_data
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
types::global_dof_index size() const
IndexSet locally_owned_range_data
bool is_compatible(const Partitioner &part) const
std::vector< std::pair< unsigned int, unsigned int > > ghost_indices_subset_data
std::map< unsigned int, IndexSet > get_requesters()
#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_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
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)
T min(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
static const unsigned int invalid_unsigned_int
#define DEAL_II_DOF_INDEX_MPI_TYPE