unsigned int
this_mpi_process(const MPI_Comm &mpi_communicator);
+ /**
+ * Return a vector of the ranks (within @p comm_large) of a subset of
+ * processes specified by @p comm_small.
+ */
+ const std::vector<unsigned int>
+ mpi_processes_within_communicator(const MPI_Comm &comm_large,
+ const MPI_Comm &comm_small);
+
/**
* Consider an unstructured communication pattern where every process in
* an MPI universe wants to send some data to a subset of the other
if (is_non_buffering_sm_supported)
{
- // TODO: move into Utilities::MPI
- const auto procs_of_sm = [](const MPI_Comm &comm,
- const MPI_Comm &comm_shared) {
- if (Utilities::MPI::job_supports_mpi() == false)
- return std::vector<unsigned int>{0};
-
- const unsigned int rank = Utilities::MPI::this_mpi_process(comm);
- const unsigned int size =
- Utilities::MPI::n_mpi_processes(comm_shared);
-
- std::vector<unsigned int> ranks(size);
- MPI_Allgather(
- &rank, 1, MPI_UNSIGNED, ranks.data(), 1, MPI_UNSIGNED, comm_shared);
-
- return ranks;
- };
-
// gather the ranks of the shared-memory domain
- const std::vector<unsigned int> sm_procs =
- procs_of_sm(task_info.communicator, communicator_sm);
+ const auto sm_procs = Utilities::MPI::mpi_processes_within_communicator(
+ task_info.communicator, communicator_sm);
// get my rank within the shared-memory domain
const unsigned int my_sm_pid = std::distance(
}
+
+ const std::vector<unsigned int>
+ mpi_processes_within_communicator(const MPI_Comm &comm_large,
+ const MPI_Comm &comm_small)
+ {
+ if (Utilities::MPI::job_supports_mpi() == false)
+ return std::vector<unsigned int>{0};
+
+ const unsigned int rank = Utilities::MPI::this_mpi_process(comm_large);
+ const unsigned int size = Utilities::MPI::n_mpi_processes(comm_small);
+
+ std::vector<unsigned int> ranks(size);
+ MPI_Allgather(
+ &rank, 1, MPI_UNSIGNED, ranks.data(), 1, MPI_UNSIGNED, comm_small);
+
+ return ranks;
+ }
+
+
+
MPI_Comm
duplicate_communicator(const MPI_Comm &mpi_communicator)
{
+ const std::vector<unsigned int>
+ mpi_processes_within_communicator(const MPI_Comm &, const MPI_Comm &)
+ {
+ return std::vector<unsigned int>{0};
+ }
+
+
+
std::vector<IndexSet>
create_ascending_partitioning(const MPI_Comm & /*comm*/,
const IndexSet::size_type local_size)