16#ifndef dealii_base_mpi_compute_index_owner_internal_h
17#define dealii_base_mpi_compute_index_owner_internal_h
36 namespace ComputeIndexOwner
47 std::pair<types::global_dof_index, types::global_dof_index>,
55 const std::map<
unsigned int,
60 const std::pair<types::global_dof_index, types::global_dof_index>
73 virtual std::vector<unsigned int>
76 std::vector<unsigned int> targets;
77 for (
const auto &rank_pair :
buffers)
78 targets.push_back(rank_pair.first);
91 &send_buffer)
override
93 send_buffer = this->
buffers.at(other_rank);
102 const unsigned int other_rank,
105 std::vector<unsigned int> &request_buffer)
override
107 (void)request_buffer;
111 for (
auto interval : buffer_recv)
137 const std::map<
unsigned int,
144 const std::pair<types::global_dof_index, types::global_dof_index>
209 std::pair<types::global_dof_index, types::global_dof_index>
248#ifdef DEAL_II_WITH_MPI
252 std::map<
unsigned int,
270 std::pair<types::global_dof_index, types::global_dof_index>
271 index_range(*interval->begin(), interval->last() + 1);
272 const unsigned int owner_last =
275 while (owner_first != owner_last)
277 Assert(index_range.first < index_range.second,
301 if (owner_first == my_rank)
308 dic_local_received += next_index - index_range.first;
313 buffers[owner_first].emplace_back(index_range.first,
316 index_range.first = next_index;
321 std::vector<MPI_Request> request;
325 owned_indices.
size())
344 for (
const auto &rank_pair : buffers)
346 request.push_back(MPI_Request());
347 const int ierr = MPI_Isend(rank_pair.second.data(),
348 rank_pair.second.size() * 2,
358 while (this->locally_owned_size != dic_local_received)
363 MPI_Probe(MPI_ANY_SOURCE, mpi_tag,
comm, &status);
368 ierr = MPI_Get_count(&status,
373 const auto other_rank = status.MPI_SOURCE;
380 buffer(number_amount / 2);
381 ierr = MPI_Recv(buffer.data(),
390 for (
auto interval : buffer)
412 dic_local_received += interval.second - interval.first;
429 std::pair<types::global_dof_index, types::global_dof_index>,
431 consensus_algo(temp,
comm);
432 consensus_algo.run();
444 if (request.size() > 0)
446 const int ierr = MPI_Waitall(request.size(),
448 MPI_STATUSES_IGNORE);
492 const unsigned int guess = 0)
501 rank_in_owned_indices);
517#ifdef DEAL_II_WITH_MPI
558 std::pair<types::global_dof_index, types::global_dof_index>,
632 std::vector<std::vector<
633 std::pair<
unsigned int,
634 std::vector<std::pair<unsigned int, unsigned int>>>>>
646 std::map<unsigned int, std::vector<types::global_dof_index>>
664 const unsigned int other_rank,
667 std::vector<unsigned int> &request_buffer)
override
669 unsigned int owner_index = 0;
670 for (
const auto &interval : buffer_recv)
671 for (
auto i = interval.first; i < interval.second; ++i)
673 const unsigned int actual_owner =
675 request_buffer.push_back(actual_owner);
686 virtual std::vector<unsigned int>
689 std::vector<unsigned int> targets;
693 unsigned int index = 0;
694 unsigned int owner_index = 0;
705 else if (targets.empty() || targets.back() != other_rank)
706 targets.push_back(other_rank);
712 for (
auto i : targets)
720 unsigned int index = 0;
748 &send_buffer)
override
753 is.
add_indices(indices_i.begin(), indices_i.end());
759 send_buffer.emplace_back(*interval->begin(),
760 interval->last() + 1);
769 const unsigned int other_rank,
770 std::vector<unsigned int> &recv_buffer)
override
781 const std::vector<unsigned int> &recv_buffer)
override
786 for (
unsigned int j = 0; j <
recv_indices[other_rank].size(); j++)
799 std::map<unsigned int, IndexSet>
803 ExcMessage(
"Must enable index range tracking in "
804 "constructor of ConsensusAlgorithmProcess"));
806 std::map<unsigned int, ::IndexSet> requested_indices;
808#ifdef DEAL_II_WITH_MPI
820 std::vector<MPI_Request> send_requests;
831 std::vector<std::vector<unsigned int>> send_data(
requesters.size());
832 for (
unsigned int i = 0; i <
requesters.size(); ++i)
841 IndexSet &my_index_set = requested_indices[j.first];
843 for (
const auto &interval : j.second)
844 my_index_set.
add_range(index_offset + interval.first,
853 send_data[i].push_back(j.first);
854 send_data[i].push_back(j.second.size());
855 for (
const auto &interval : j.second)
857 send_data[i].push_back(interval.first);
858 send_data[i].push_back(interval.second);
861 send_requests.push_back(MPI_Request());
863 MPI_Isend(send_data[i].data(),
869 &send_requests.back());
881 MPI_Probe(MPI_ANY_SOURCE, mpi_tag,
comm, &status);
886 ierr = MPI_Get_count(&status, MPI_UNSIGNED, &number_amount);
891 std::vector<std::pair<unsigned int, unsigned int>> buffer(
893 ierr = MPI_Recv(buffer.data(),
906 unsigned int offset = 0;
907 while (offset < buffer.size())
913 for (
unsigned int i = offset + 1;
914 i < offset + buffer[offset].second + 1;
916 my_index_set.
add_range(index_offset + buffer[i].first,
917 index_offset + buffer[i].second);
923 requested_indices[buffer[offset].first];
924 if (index_set.
size() == 0)
928 offset += buffer[offset].second + 1;
933 if (send_requests.size() > 0)
935 const auto ierr = MPI_Waitall(send_requests.size(),
936 send_requests.data(),
937 MPI_STATUSES_IGNORE);
943 for (
const auto &it : requested_indices)
949 "The indices requested from the current "
950 "MPI rank should be locally owned here!"));
956 return requested_indices;
975 unsigned int & owner_index,
976 const unsigned int rank_of_request)
982 const unsigned int rank_of_owner =
990 std::vector<std::pair<unsigned int, unsigned int>>());
992 requesters[owner_index].back().second.back().second !=
994 requesters[owner_index].back().second.emplace_back(
998 ++
requesters[owner_index].back().second.back().second;
IntervalIterator end_intervals() const
size_type n_elements() const
void set_size(const size_type size)
IntervalIterator begin_intervals() const
void subtract_set(const IndexSet &other)
void add_range(const size_type begin, const size_type end)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
const unsigned int my_rank
virtual void prepare_buffer_for_answer(const unsigned int other_rank, std::vector< unsigned int > &recv_buffer) override
std::vector< std::vector< std::pair< unsigned int, std::vector< std::pair< unsigned int, unsigned int > > > > > requesters
std::vector< unsigned int > & owning_ranks
void append_index_origin(const types::global_dof_index index, unsigned int &owner_index, const unsigned int rank_of_request)
virtual std::vector< unsigned int > compute_targets() override
const IndexSet & indices_to_look_up
virtual void read_answer(const unsigned int other_rank, const std::vector< unsigned int > &recv_buffer) override
const unsigned int n_procs
const IndexSet & owned_indices
std::map< unsigned int, IndexSet > get_requesters()
std::map< unsigned int, std::vector< unsigned int > > recv_indices
virtual void answer_request(const unsigned int other_rank, const std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &buffer_recv, std::vector< unsigned int > &request_buffer) override
virtual void create_request(const unsigned int other_rank, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &send_buffer) override
std::map< unsigned int, std::vector< types::global_dof_index > > indices_to_look_up_by_dict_rank
ConsensusAlgorithmsPayload(const IndexSet &owned_indices, const IndexSet &indices_to_look_up, const MPI_Comm &comm, std::vector< unsigned int > &owning_ranks, const bool track_index_requests=false)
const bool track_index_requests
DictionaryPayLoad(const std::map< unsigned int, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > &buffers, std::vector< unsigned int > &actually_owning_ranks, const std::pair< types::global_dof_index, types::global_dof_index > &local_range, std::vector< unsigned int > &actually_owning_rank_list)
virtual std::vector< unsigned int > compute_targets() override
const std::map< unsigned int, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > & buffers
std::vector< unsigned int > & actually_owning_rank_list
std::vector< unsigned int > & actually_owning_ranks
virtual void create_request(const unsigned int other_rank, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &send_buffer) override
const std::pair< types::global_dof_index, types::global_dof_index > & local_range
virtual void answer_request(const unsigned int other_rank, const std::vector< std::pair< types::global_dof_index, types::global_dof_index > > &buffer_recv, std::vector< unsigned int > &request_buffer) override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int this_mpi_process(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)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
*braid_SplitCommworld & comm
std::vector< unsigned int > actually_owning_rank_list
void reinit(const IndexSet &owned_indices, const MPI_Comm &comm)
unsigned int dof_to_dict_rank(const types::global_dof_index i)
types::global_dof_index get_index_offset(const unsigned int rank)
types::global_dof_index locally_owned_size
types::global_dof_index dofs_per_process
unsigned int n_dict_procs_in_owned_indices
unsigned int stride_small_size
std::pair< types::global_dof_index, types::global_dof_index > local_range
std::vector< unsigned int > actually_owning_ranks
void partition(const IndexSet &owned_indices, const MPI_Comm &comm)
static const unsigned int range_minimum_grain_size
unsigned int get_owning_rank_index(const unsigned int rank_in_owned_indices, const unsigned int guess=0)
types::global_dof_index size
#define DEAL_II_DOF_INDEX_MPI_TYPE