Reference documentation for deal.II version 9.2.0
\(\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\}}\)
mpi_compute_index_owner_internal.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2019 - 2020 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 
16 #ifndef dealii_base_mpi_compute_index_owner_internal_h
17 #define dealii_base_mpi_compute_index_owner_internal_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/mpi.h>
23 
25 
26 namespace Utilities
27 {
28  namespace MPI
29  {
30  namespace internal
31  {
36  namespace ComputeIndexOwner
37  {
47  std::pair<types::global_dof_index, types::global_dof_index>,
48  unsigned int>
49  {
50  public:
55  const std::map<unsigned int,
56  std::vector<std::pair<types::global_dof_index,
58  & buffers,
59  std::vector<unsigned int> &actually_owning_ranks,
60  const std::pair<types::global_dof_index, types::global_dof_index>
61  & local_range,
62  std::vector<unsigned int> &actually_owning_rank_list)
63  : buffers(buffers)
67  {}
68 
73  virtual std::vector<unsigned int>
74  compute_targets() override
75  {
76  std::vector<unsigned int> targets;
77  for (const auto &rank_pair : buffers)
78  targets.push_back(rank_pair.first);
79 
80  return targets;
81  }
82 
87  virtual void
88  create_request(const unsigned int other_rank,
89  std::vector<std::pair<types::global_dof_index,
91  &send_buffer) override
92  {
93  send_buffer = this->buffers.at(other_rank);
94  }
95 
100  virtual void
102  const unsigned int other_rank,
103  const std::vector<std::pair<types::global_dof_index,
104  types::global_dof_index>> &buffer_recv,
105  std::vector<unsigned int> &request_buffer) override
106  {
107  (void)request_buffer; // not needed
108 
109 
110  // process message: loop over all intervals
111  for (auto interval : buffer_recv)
112  {
113 #ifdef DEBUG
114  for (types::global_dof_index i = interval.first;
115  i < interval.second;
116  i++)
119  ExcInternalError());
120  Assert(interval.first >= local_range.first &&
121  interval.first < local_range.second,
122  ExcInternalError());
123  Assert(interval.second > local_range.first &&
124  interval.second <= local_range.second,
125  ExcInternalError());
126 #endif
127  std::fill(actually_owning_ranks.data() + interval.first -
128  local_range.first,
129  actually_owning_ranks.data() + interval.second -
130  local_range.first,
131  other_rank);
132  }
133  actually_owning_rank_list.push_back(other_rank);
134  }
135 
136  private:
137  const std::map<unsigned int,
138  std::vector<std::pair<types::global_dof_index,
141 
142  std::vector<unsigned int> &actually_owning_ranks;
143 
144  const std::pair<types::global_dof_index, types::global_dof_index>
146 
147  std::vector<unsigned int> &actually_owning_rank_list;
148  };
149 
150 
151 
157  struct Dictionary
158  {
162  static const unsigned int range_minimum_grain_size = 4096;
163 
170  std::vector<unsigned int> actually_owning_ranks;
171 
176  std::vector<unsigned int> actually_owning_rank_list;
177 
185 
191  std::pair<types::global_dof_index, types::global_dof_index>
193 
200 
205 
211 
217  unsigned int stride_small_size;
218 
224  void
225  reinit(const IndexSet &owned_indices, const MPI_Comm &comm)
226  {
227  // 1) set up the partition
228  this->partition(owned_indices, comm);
229 
230 #ifdef DEAL_II_WITH_MPI
231  unsigned int my_rank = this_mpi_process(comm);
232 
233  types::global_dof_index dic_local_received = 0;
234  std::map<unsigned int,
235  std::vector<std::pair<types::global_dof_index,
237  buffers;
238 
239  std::fill(actually_owning_ranks.begin(),
240  actually_owning_ranks.end(),
242 
243  // 2) collect relevant processes and process local dict entries
244  for (auto interval = owned_indices.begin_intervals();
245  interval != owned_indices.end_intervals();
246  ++interval)
247  {
248  // Due to the granularity of the dictionary, the interval
249  // might be split into several ranges of processor owner
250  // ranks. Here, we process the interval by breaking into
251  // smaller pieces in terms of the dictionary number.
252  std::pair<types::global_dof_index, types::global_dof_index>
253  index_range(*interval->begin(), interval->last() + 1);
254  const unsigned int owner_last =
255  dof_to_dict_rank(interval->last());
256  unsigned int owner_first = numbers::invalid_unsigned_int;
257  while (owner_first != owner_last)
258  {
259  Assert(index_range.first < index_range.second,
260  ExcInternalError());
261 
262  owner_first = dof_to_dict_rank(index_range.first);
263 
264  // this explicitly picks up the formula of
265  // dof_to_dict_rank, so the two places must be in sync
266  types::global_dof_index next_index =
267  std::min(get_index_offset(owner_first + 1),
268  index_range.second);
269 
270  Assert(next_index > index_range.first, ExcInternalError());
271 
272 # ifdef DEBUG
273  // make sure that the owner is the same on the current
274  // interval
275  for (types::global_dof_index i = index_range.first + 1;
276  i < next_index;
277  ++i)
278  AssertDimension(owner_first, dof_to_dict_rank(i));
279 # endif
280 
281  // add the interval, either to the local range or into a
282  // buffer to be sent to another processor
283  if (owner_first == my_rank)
284  {
285  std::fill(actually_owning_ranks.data() +
286  index_range.first - local_range.first,
287  actually_owning_ranks.data() + next_index -
288  local_range.first,
289  my_rank);
290  dic_local_received += next_index - index_range.first;
291  if (actually_owning_rank_list.empty())
292  actually_owning_rank_list.push_back(my_rank);
293  }
294  else
295  buffers[owner_first].emplace_back(index_range.first,
296  next_index);
297 
298  index_range.first = next_index;
299  }
300  }
301 
302  n_dict_procs_in_owned_indices = buffers.size();
303  std::vector<MPI_Request> request;
304 
305  // Check if index set space is partitioned globally without gaps.
306  if (Utilities::MPI::sum(owned_indices.n_elements(), comm) ==
307  owned_indices.size())
308  {
309  // no gaps: setup is simple! Processes send their locally owned
310  // indices to the dictionary. The dictionary stores the sending
311  // rank for each index. The dictionary knows exactly
312  // when it is set up when all indices it is responsible for
313  // have been processed.
314 
315  request.reserve(n_dict_procs_in_owned_indices);
316 
317  // protect the following communication steps using a mutex:
318  static CollectiveMutex mutex;
319  CollectiveMutex::ScopedLock lock(mutex, comm);
320 
321  const int mpi_tag =
323 
324 
325  // 3) send messages with local dofs to the right dict process
326  for (const auto &rank_pair : buffers)
327  {
328  request.push_back(MPI_Request());
329  const int ierr = MPI_Isend(rank_pair.second.data(),
330  rank_pair.second.size() * 2,
332  rank_pair.first,
333  mpi_tag,
334  comm,
335  &request.back());
336  AssertThrowMPI(ierr);
337  }
338 
339  // 4) receive messages until all dofs in dict are processed
340  while (this->local_size != dic_local_received)
341  {
342  // wait for an incoming message
343  MPI_Status status;
344  auto ierr =
345  MPI_Probe(MPI_ANY_SOURCE, mpi_tag, comm, &status);
346  AssertThrowMPI(ierr);
347 
348  // retrieve size of incoming message
349  int number_amount;
350  ierr = MPI_Get_count(&status,
352  &number_amount);
353  AssertThrowMPI(ierr);
354 
355  const auto other_rank = status.MPI_SOURCE;
356  actually_owning_rank_list.push_back(other_rank);
357 
358  // receive message
359  Assert(number_amount % 2 == 0, ExcInternalError());
360  std::vector<std::pair<types::global_dof_index,
362  buffer(number_amount / 2);
363  ierr = MPI_Recv(buffer.data(),
364  number_amount,
366  status.MPI_SOURCE,
367  status.MPI_TAG,
368  comm,
369  MPI_STATUS_IGNORE);
370  AssertThrowMPI(ierr);
371  // process message: loop over all intervals
372  for (auto interval : buffer)
373  {
374 # ifdef DEBUG
375  for (types::global_dof_index i = interval.first;
376  i < interval.second;
377  i++)
380  ExcInternalError());
381  Assert(interval.first >= local_range.first &&
382  interval.first < local_range.second,
383  ExcInternalError());
384  Assert(interval.second > local_range.first &&
385  interval.second <= local_range.second,
386  ExcInternalError());
387 # endif
388 
389  std::fill(actually_owning_ranks.data() +
390  interval.first - local_range.first,
391  actually_owning_ranks.data() +
392  interval.second - local_range.first,
393  other_rank);
394  dic_local_received += interval.second - interval.first;
395  }
396  }
397  }
398  else
399  {
400  // with gap: use a ConsensusAlgorithm to determine when all
401  // dictionaries have been set up.
402 
403  // 3/4) use a ConsensusAlgorithm to send messages with local
404  // dofs to the right dict process
405  DictionaryPayLoad temp(buffers,
407  local_range,
409 
411  std::pair<types::global_dof_index, types::global_dof_index>,
412  unsigned int>
413  consensus_algo(temp, comm);
414  consensus_algo.run();
415  }
416 
417  std::sort(actually_owning_rank_list.begin(),
419 
420  for (unsigned int i = 1; i < actually_owning_rank_list.size(); ++i)
423  ExcInternalError());
424 
425  // 5) make sure that all messages have been sent
426  if (request.size() > 0)
427  {
428  const int ierr = MPI_Waitall(request.size(),
429  request.data(),
430  MPI_STATUSES_IGNORE);
431  AssertThrowMPI(ierr);
432  }
433 
434 #else
435  (void)owned_indices;
436  (void)comm;
437 #endif
438  }
439 
445  unsigned int
447  {
448  // note: this formula is also explicitly used in
449  // get_index_offset(), so keep the two in sync
450  return (i / dofs_per_process) * stride_small_size;
451  }
452 
458  get_index_offset(const unsigned int rank)
459  {
460  return std::min(dofs_per_process *
461  static_cast<types::global_dof_index>(
462  (rank + stride_small_size - 1) /
464  size);
465  }
466 
472  unsigned int
473  get_owning_rank_index(const unsigned int rank_in_owned_indices,
474  const unsigned int guess = 0)
475  {
477  if (actually_owning_rank_list[guess] == rank_in_owned_indices)
478  return guess;
479  else
480  {
483  rank_in_owned_indices);
484  Assert(it != actually_owning_rank_list.end(),
485  ExcInternalError());
486  Assert(*it == rank_in_owned_indices, ExcInternalError());
487  return it - actually_owning_rank_list.begin();
488  }
489  }
490 
491  private:
496  void
497  partition(const IndexSet &owned_indices, const MPI_Comm &comm)
498  {
499 #ifdef DEAL_II_WITH_MPI
500  const unsigned int n_procs = n_mpi_processes(comm);
501  const unsigned int my_rank = this_mpi_process(comm);
502 
503  size = owned_indices.size();
504 
505  Assert(size > 0, ExcNotImplemented());
506 
507  dofs_per_process = (size + n_procs - 1) / n_procs;
509  {
512  }
513  else
514  stride_small_size = 1;
515  local_range.first = get_index_offset(my_rank);
516  local_range.second = get_index_offset(my_rank + 1);
517 
518  local_size = local_range.second - local_range.first;
519 
523 #else
524  (void)owned_indices;
525  (void)comm;
526 #endif
527  }
528  };
529 
530 
531 
540  std::pair<types::global_dof_index, types::global_dof_index>,
541  unsigned int>
542  {
543  public:
549  const MPI_Comm &comm,
550  std::vector<unsigned int> &owning_ranks,
551  const bool track_index_requests = false)
554  , comm(comm)
559  {
562  }
563 
568 
574 
578  const MPI_Comm comm;
579 
583  const unsigned int my_rank;
584 
589  const unsigned int n_procs;
590 
597 
603  std::vector<unsigned int> &owning_ranks;
604 
614  std::vector<std::vector<
615  std::pair<unsigned int,
616  std::vector<std::pair<unsigned int, unsigned int>>>>>
618 
623 
628  std::map<unsigned int, std::vector<types::global_dof_index>>
630 
635  std::map<unsigned int, std::vector<unsigned int>> recv_indices;
636 
644  virtual void
646  const unsigned int other_rank,
647  const std::vector<std::pair<types::global_dof_index,
648  types::global_dof_index>> &buffer_recv,
649  std::vector<unsigned int> &request_buffer) override
650  {
651  unsigned int owner_index = 0;
652  for (const auto &interval : buffer_recv)
653  for (auto i = interval.first; i < interval.second; ++i)
654  {
655  const unsigned int actual_owner =
657  request_buffer.push_back(actual_owner);
658 
660  append_index_origin(i, owner_index, other_rank);
661  }
662  }
663 
668  virtual std::vector<unsigned int>
669  compute_targets() override
670  {
671  std::vector<unsigned int> targets;
672 
673  // 1) collect relevant processes and process local dict entries
674  {
675  unsigned int index = 0;
676  unsigned int owner_index = 0;
677  for (auto i : indices_to_look_up)
678  {
679  unsigned int other_rank = dict.dof_to_dict_rank(i);
680  if (other_rank == my_rank)
681  {
682  owning_ranks[index] =
685  append_index_origin(i, owner_index, my_rank);
686  }
687  else if (targets.empty() || targets.back() != other_rank)
688  targets.push_back(other_rank);
689  index++;
690  }
691  }
692 
693 
694  for (auto i : targets)
695  {
696  recv_indices[i] = {};
698  }
699 
700  // 3) collect indices for each process
701  {
702  unsigned int index = 0;
703  for (auto i : indices_to_look_up)
704  {
705  unsigned int other_rank = dict.dof_to_dict_rank(i);
706  if (other_rank != my_rank)
707  {
708  recv_indices[other_rank].push_back(index);
709  indices_to_look_up_by_dict_rank[other_rank].push_back(i);
710  }
711  index++;
712  }
713  }
714 
715  Assert(targets.size() == recv_indices.size() &&
716  targets.size() == indices_to_look_up_by_dict_rank.size(),
717  ExcMessage("Size does not match!"));
718 
719  return targets;
720  }
721 
726  virtual void
727  create_request(const unsigned int other_rank,
728  std::vector<std::pair<types::global_dof_index,
730  &send_buffer) override
731  {
732  // create index set and compress data to be sent
733  auto & indices_i = indices_to_look_up_by_dict_rank[other_rank];
734  IndexSet is(dict.size);
735  is.add_indices(indices_i.begin(), indices_i.end());
736  is.compress();
737 
738  for (auto interval = is.begin_intervals();
739  interval != is.end_intervals();
740  ++interval)
741  send_buffer.emplace_back(*interval->begin(),
742  interval->last() + 1);
743  }
744 
749  virtual void
751  const unsigned int other_rank,
752  std::vector<unsigned int> &recv_buffer) override
753  {
754  recv_buffer.resize(recv_indices[other_rank].size());
755  }
756 
761  virtual void
762  read_answer(const unsigned int other_rank,
763  const std::vector<unsigned int> &recv_buffer) override
764  {
765  Assert(recv_indices[other_rank].size() == recv_buffer.size(),
766  ExcMessage("Sizes do not match!"));
767 
768  for (unsigned int j = 0; j < recv_indices[other_rank].size(); j++)
769  owning_ranks[recv_indices[other_rank][j]] = recv_buffer[j];
770  }
771 
781  std::map<unsigned int, IndexSet>
783  {
785  ExcMessage("Must enable index range tracking in "
786  "constructor of ConsensusAlgorithmProcess"));
787 
788  std::map<unsigned int, ::IndexSet> requested_indices;
789 
790 #ifdef DEAL_II_WITH_MPI
791 
792  static CollectiveMutex mutex;
793  CollectiveMutex::ScopedLock lock(mutex, comm);
794 
795  const int mpi_tag = Utilities::MPI::internal::Tags::
797 
798  // reserve enough slots for the requests ahead; depending on
799  // whether the owning rank is one of the requesters or not, we
800  // might have one less requests to execute, so fill the requests
801  // on demand.
802  std::vector<MPI_Request> send_requests;
803  send_requests.reserve(requesters.size());
804 
805  // We use an integer vector for the data exchange. Since we send
806  // data associated to intervals with different requesters, we will
807  // need to send (a) the MPI rank of the requester, (b) the number
808  // of intervals directed to this requester, and (c) a list of
809  // intervals, i.e., two integers per interval. The number of items
810  // sent in total can be deduced both via the MPI status message at
811  // the receiver site as well as be counting the buckets from
812  // different requesters.
813  std::vector<std::vector<unsigned int>> send_data(requesters.size());
814  for (unsigned int i = 0; i < requesters.size(); ++i)
815  {
816  // special code for our own indices
818  {
819  for (const auto &j : requesters[i])
820  {
821  const types::global_dof_index index_offset =
823  IndexSet &my_index_set = requested_indices[j.first];
824  my_index_set.set_size(owned_indices.size());
825  for (const auto &interval : j.second)
826  my_index_set.add_range(index_offset + interval.first,
827  index_offset +
828  interval.second);
829  }
830  }
831  else
832  {
833  for (const auto &j : requesters[i])
834  {
835  send_data[i].push_back(j.first);
836  send_data[i].push_back(j.second.size());
837  for (const auto &interval : j.second)
838  {
839  send_data[i].push_back(interval.first);
840  send_data[i].push_back(interval.second);
841  }
842  }
843  send_requests.push_back(MPI_Request());
844  const int ierr =
845  MPI_Isend(send_data[i].data(),
846  send_data[i].size(),
847  MPI_UNSIGNED,
849  mpi_tag,
850  comm,
851  &send_requests.back());
852  AssertThrowMPI(ierr);
853  }
854  }
855 
856  // receive the data
857  for (unsigned int c = 0; c < dict.n_dict_procs_in_owned_indices;
858  ++c)
859  {
860  // wait for an incoming message
861  MPI_Status status;
862  unsigned int ierr =
863  MPI_Probe(MPI_ANY_SOURCE, mpi_tag, comm, &status);
864  AssertThrowMPI(ierr);
865 
866  // retrieve size of incoming message
867  int number_amount;
868  ierr = MPI_Get_count(&status, MPI_UNSIGNED, &number_amount);
869  AssertThrowMPI(ierr);
870 
871  // receive message
872  Assert(number_amount % 2 == 0, ExcInternalError());
873  std::vector<std::pair<unsigned int, unsigned int>> buffer(
874  number_amount / 2);
875  ierr = MPI_Recv(buffer.data(),
876  number_amount,
877  MPI_UNSIGNED,
878  status.MPI_SOURCE,
879  status.MPI_TAG,
880  comm,
881  &status);
882  AssertThrowMPI(ierr);
883 
884  // unpack the message and translate the dictionary-local
885  // indices coming via MPI to the global index range
886  const types::global_dof_index index_offset =
887  dict.get_index_offset(status.MPI_SOURCE);
888  unsigned int offset = 0;
889  while (offset < buffer.size())
890  {
891  AssertIndexRange(offset + buffer[offset].second,
892  buffer.size());
893 
894  IndexSet my_index_set(owned_indices.size());
895  for (unsigned int i = offset + 1;
896  i < offset + buffer[offset].second + 1;
897  ++i)
898  my_index_set.add_range(index_offset + buffer[i].first,
899  index_offset + buffer[i].second);
900 
901  // the underlying index set is able to merge ranges coming
902  // from different ranks due to the partitioning in the
903  // dictionary
904  IndexSet &index_set =
905  requested_indices[buffer[offset].first];
906  if (index_set.size() == 0)
907  index_set.set_size(owned_indices.size());
908  index_set.add_indices(my_index_set);
909 
910  offset += buffer[offset].second + 1;
911  }
912  AssertDimension(offset, buffer.size());
913  }
914 
915  if (send_requests.size() > 0)
916  {
917  const auto ierr = MPI_Waitall(send_requests.size(),
918  send_requests.data(),
919  MPI_STATUSES_IGNORE);
920  AssertThrowMPI(ierr);
921  }
922 
923 
924 # ifdef DEBUG
925  for (const auto &it : requested_indices)
926  {
927  IndexSet copy_set = it.second;
928  copy_set.subtract_set(owned_indices);
929  Assert(copy_set.n_elements() == 0,
931  "The indices requested from the current "
932  "MPI rank should be locally owned here!"));
933  }
934 # endif
935 
936 #endif // DEAL_II_WITH_MPI
937 
938  return requested_indices;
939  }
940 
941  private:
955  void
957  unsigned int & owner_index,
958  const unsigned int rank_of_request)
959  {
960  // remember who requested which index. We want to use an
961  // std::vector with simple addressing, via a good guess from the
962  // preceding index, rather than std::map, because this is an inner
963  // loop and it avoids the map lookup in every iteration
964  const unsigned int rank_of_owner =
966  owner_index =
967  dict.get_owning_rank_index(rank_of_owner, owner_index);
968  if (requesters[owner_index].empty() ||
969  requesters[owner_index].back().first != rank_of_request)
970  requesters[owner_index].emplace_back(
971  rank_of_request,
972  std::vector<std::pair<unsigned int, unsigned int>>());
973  if (requesters[owner_index].back().second.empty() ||
974  requesters[owner_index].back().second.back().second !=
975  index - dict.local_range.first)
976  requesters[owner_index].back().second.emplace_back(
977  index - dict.local_range.first,
978  index - dict.local_range.first + 1);
979  else
980  ++requesters[owner_index].back().second.back().second;
981  }
982  };
983 
984  } // namespace ComputeIndexOwner
985  } // namespace internal
986  } // namespace MPI
987 } // namespace Utilities
988 
990 
991 #endif
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::compute_targets
virtual std::vector< unsigned int > compute_targets() override
Definition: mpi_compute_index_owner_internal.h:74
IndexSet::compress
void compress() const
Definition: index_set.h:1643
IndexSet
Definition: index_set.h:74
IndexSet::subtract_set
void subtract_set(const IndexSet &other)
Definition: index_set.cc:258
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::create_request
virtual void create_request(const unsigned int other_rank, std::vector< std::pair< types::global_dof_index, types::global_dof_index >> &send_buffer) override
Definition: mpi_compute_index_owner_internal.h:727
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::comm
const MPI_Comm comm
Definition: mpi_compute_index_owner_internal.h:578
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::actually_owning_ranks
std::vector< unsigned int > actually_owning_ranks
Definition: mpi_compute_index_owner_internal.h:170
DEAL_II_DOF_INDEX_MPI_TYPE
#define DEAL_II_DOF_INDEX_MPI_TYPE
Definition: types.h:86
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::requesters
std::vector< std::vector< std::pair< unsigned int, std::vector< std::pair< unsigned int, unsigned int > > > > > requesters
Definition: mpi_compute_index_owner_internal.h:617
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::get_requesters
std::map< unsigned int, IndexSet > get_requesters()
Definition: mpi_compute_index_owner_internal.h:782
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
Definition: mpi_compute_index_owner_internal.h:538
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::create_request
virtual void create_request(const unsigned int other_rank, std::vector< std::pair< types::global_dof_index, types::global_dof_index >> &send_buffer) override
Definition: mpi_compute_index_owner_internal.h:88
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::partition
void partition(const IndexSet &owned_indices, const MPI_Comm &comm)
Definition: mpi_compute_index_owner_internal.h:497
IndexSet::size
size_type size() const
Definition: index_set.h:1635
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::actually_owning_rank_list
std::vector< unsigned int > & actually_owning_rank_list
Definition: mpi_compute_index_owner_internal.h:147
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::owned_indices
const IndexSet & owned_indices
Definition: mpi_compute_index_owner_internal.h:567
MPI_Comm
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::ConsensusAlgorithmsPayload
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)
Definition: mpi_compute_index_owner_internal.h:547
second
Point< 2 > second
Definition: grid_out.cc:4353
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::stride_small_size
unsigned int stride_small_size
Definition: mpi_compute_index_owner_internal.h:217
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::compute_targets
virtual std::vector< unsigned int > compute_targets() override
Definition: mpi_compute_index_owner_internal.h:669
mpi_consensus_algorithms.h
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::recv_indices
std::map< unsigned int, std::vector< unsigned int > > recv_indices
Definition: mpi_compute_index_owner_internal.h:635
IndexSet::begin_intervals
IntervalIterator begin_intervals() const
Definition: index_set.h:1590
Utilities::MPI::ConsensusAlgorithms::Selector
Definition: mpi_consensus_algorithms.h:464
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::size
types::global_dof_index size
Definition: mpi_compute_index_owner_internal.h:204
IndexSet::n_elements
size_type n_elements() const
Definition: index_set.h:1833
IndexSet::set_size
void set_size(const size_type size)
Definition: index_set.h:1623
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::indices_to_look_up
const IndexSet & indices_to_look_up
Definition: mpi_compute_index_owner_internal.h:573
IndexSet::end_intervals
IntervalIterator end_intervals() const
Definition: index_set.h:1602
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::dof_to_dict_rank
unsigned int dof_to_dict_rank(const types::global_dof_index i)
Definition: mpi_compute_index_owner_internal.h:446
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::dofs_per_process
types::global_dof_index dofs_per_process
Definition: mpi_compute_index_owner_internal.h:184
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::actually_owning_ranks
std::vector< unsigned int > & actually_owning_ranks
Definition: mpi_compute_index_owner_internal.h:142
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::reinit
void reinit(const IndexSet &owned_indices, const MPI_Comm &comm)
Definition: mpi_compute_index_owner_internal.h:225
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::dict
Dictionary dict
Definition: mpi_compute_index_owner_internal.h:622
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::local_range
std::pair< types::global_dof_index, types::global_dof_index > local_range
Definition: mpi_compute_index_owner_internal.h:192
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::indices_to_look_up_by_dict_rank
std::map< unsigned int, std::vector< types::global_dof_index > > indices_to_look_up_by_dict_rank
Definition: mpi_compute_index_owner_internal.h:629
Utilities::MPI::internal::ComputeIndexOwner::Dictionary
Definition: mpi_compute_index_owner_internal.h:157
mpi.h
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::my_rank
const unsigned int my_rank
Definition: mpi_compute_index_owner_internal.h:583
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad
Definition: mpi_compute_index_owner_internal.h:45
Utilities::lower_bound
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1102
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::track_index_requests
const bool track_index_requests
Definition: mpi_compute_index_owner_internal.h:596
Utilities::MPI::CollectiveMutex
Definition: mpi.h:322
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::buffers
const std::map< unsigned int, std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > & buffers
Definition: mpi_compute_index_owner_internal.h:140
Utilities::MPI::internal::Tags::dictionary_reinit
@ dictionary_reinit
Dictionary::reinit()
Definition: mpi_tags.h:78
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::answer_request
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
Definition: mpi_compute_index_owner_internal.h:101
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::answer_request
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
Definition: mpi_compute_index_owner_internal.h:645
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::n_procs
const unsigned int n_procs
Definition: mpi_compute_index_owner_internal.h:589
unsigned int
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::actually_owning_rank_list
std::vector< unsigned int > actually_owning_rank_list
Definition: mpi_compute_index_owner_internal.h:176
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::owning_ranks
std::vector< unsigned int > & owning_ranks
Definition: mpi_compute_index_owner_internal.h:603
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::read_answer
virtual void read_answer(const unsigned int other_rank, const std::vector< unsigned int > &recv_buffer) override
Definition: mpi_compute_index_owner_internal.h:762
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::DictionaryPayLoad
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)
Definition: mpi_compute_index_owner_internal.h:54
Utilities::MPI::this_mpi_process
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
Utilities::MPI::internal::ComputeIndexOwner::DictionaryPayLoad::local_range
const std::pair< types::global_dof_index, types::global_dof_index > & local_range
Definition: mpi_compute_index_owner_internal.h:145
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::range_minimum_grain_size
static const unsigned int range_minimum_grain_size
Definition: mpi_compute_index_owner_internal.h:162
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
Utilities::MPI::internal::Tags::consensus_algorithm_payload_get_requesters
@ consensus_algorithm_payload_get_requesters
ConsensusAlgorithms::Payload::get_requesters()
Definition: mpi_tags.h:81
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::append_index_origin
void append_index_origin(const types::global_dof_index index, unsigned int &owner_index, const unsigned int rank_of_request)
Definition: mpi_compute_index_owner_internal.h:956
config.h
numbers::invalid_subdomain_id
const types::subdomain_id invalid_subdomain_id
Definition: types.h:285
internal
Definition: aligned_vector.h:369
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::n_dict_procs_in_owned_indices
unsigned int n_dict_procs_in_owned_indices
Definition: mpi_compute_index_owner_internal.h:210
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
IndexSet::add_indices
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1704
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::get_index_offset
types::global_dof_index get_index_offset(const unsigned int rank)
Definition: mpi_compute_index_owner_internal.h:458
first
Point< 2 > first
Definition: grid_out.cc:4352
Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload::prepare_buffer_for_answer
virtual void prepare_buffer_for_answer(const unsigned int other_rank, std::vector< unsigned int > &recv_buffer) override
Definition: mpi_compute_index_owner_internal.h:750
MPI_Request
Utilities
Definition: cuda.h:31
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::local_size
types::global_dof_index local_size
Definition: mpi_compute_index_owner_internal.h:199
Utilities::MPI::internal::ComputeIndexOwner::Dictionary::get_owning_rank_index
unsigned int get_owning_rank_index(const unsigned int rank_in_owned_indices, const unsigned int guess=0)
Definition: mpi_compute_index_owner_internal.h:473
IndexSet::add_range
void add_range(const size_type begin, const size_type end)
Definition: index_set.h:1674
Utilities::MPI::CollectiveMutex::ScopedLock
Definition: mpi.h:330
int
Utilities::MPI::ConsensusAlgorithms::Process
Definition: mpi_consensus_algorithms.h:64