template <typename T>
- std::vector<T>
- compute_set_union(const std::vector<T> &vec, const MPI_Comm &comm)
+ T
+ all_reduce(const T & vec,
+ const MPI_Comm & comm,
+ const std::function<T(const T &, const T &)> &combiner)
{
#ifdef DEAL_II_WITH_MPI
- // 1) collect vector entries and create union
- std::vector<T> result = vec;
-
- const unsigned int rank = this_mpi_process(comm);
- const unsigned int nproc = n_mpi_processes(comm);
-
- for (unsigned int stride = 1; stride < nproc; stride *= 2)
+ if (job_supports_mpi())
{
- const unsigned int rank_recv = (2 * stride) * (rank / (2 * stride));
- const unsigned int rank_send = rank_recv + stride;
+ // 1) perform custom reduction
+ T result = vec;
- if (rank_send >= nproc) // nothing to do
- continue;
+ const unsigned int rank = this_mpi_process(comm);
+ const unsigned int nproc = n_mpi_processes(comm);
- if (rank_recv == rank) // process receives data
+ for (unsigned int stride = 1; stride < nproc; stride *= 2)
{
- MPI_Status status;
- const auto ierr_1 = MPI_Probe(rank_send,
- internal::Tags::compute_union,
- comm,
- &status);
- AssertThrowMPI(ierr_1);
-
- int amount;
- const auto ierr_2 =
- MPI_Get_count(&status,
- internal::mpi_type_id(result.data()),
- &amount);
- AssertThrowMPI(ierr_2);
-
- const unsigned int old_size = result.size();
- result.resize(old_size + amount);
-
- const auto ierr_3 = MPI_Recv(result.data() + old_size,
- amount,
- internal::mpi_type_id(result.data()),
- status.MPI_SOURCE,
- status.MPI_TAG,
- comm,
- MPI_STATUS_IGNORE);
- AssertThrowMPI(ierr_3);
-
- std::sort(result.begin(), result.end());
- result.erase(std::unique(result.begin(), result.end()),
- result.end());
+ const unsigned int rank_recv =
+ (2 * stride) * (rank / (2 * stride));
+ const unsigned int rank_send = rank_recv + stride;
+
+ if (rank_send >= nproc) // nothing to do
+ continue;
+
+ if (rank_recv == rank) // process receives data
+ {
+ MPI_Status status;
+ const auto ierr_1 = MPI_Probe(rank_send,
+ internal::Tags::compute_union,
+ comm,
+ &status);
+ AssertThrowMPI(ierr_1);
+
+ int amount;
+ const auto ierr_2 = MPI_Get_count(&status, MPI_CHAR, &amount);
+ AssertThrowMPI(ierr_2);
+
+ std::vector<char> temp(amount);
+
+ const auto ierr_3 = MPI_Recv(temp.data(),
+ amount,
+ MPI_CHAR,
+ status.MPI_SOURCE,
+ status.MPI_TAG,
+ comm,
+ MPI_STATUS_IGNORE);
+ AssertThrowMPI(ierr_3);
+
+ result = combiner(result, Utilities::unpack<T>(temp, false));
+ }
+ else if (rank_send == rank) // process sends data
+ {
+ const std::vector<char> temp = Utilities::pack(result, false);
+
+ const auto ierr = MPI_Send(temp.data(),
+ temp.size(),
+ MPI_CHAR,
+ rank_recv,
+ internal::Tags::compute_union,
+ comm);
+ AssertThrowMPI(ierr);
+ }
}
- else if (rank_send == rank) // process sends data
- {
- const auto ierr = MPI_Send(result.data(),
- result.size(),
- internal::mpi_type_id(result.data()),
- rank_recv,
- internal::Tags::compute_union,
- comm);
- AssertThrowMPI(ierr);
- }
- }
- // 2) broadcast result
- int size = result.size();
- MPI_Bcast(&size, 1, MPI_INT, 0, comm);
- result.resize(size);
- MPI_Bcast(
- result.data(), size, internal::mpi_type_id(vec.data()), 0, comm);
+ // 2) broadcast result
+ std::vector<char> temp = Utilities::pack(result, false);
+ unsigned int size = temp.size();
+ MPI_Bcast(&size, 1, MPI_UNSIGNED, 0, comm);
+ temp.resize(size);
+ MPI_Bcast(temp.data(), size, MPI_CHAR, 0, comm);
- return result;
-#else
+ return Utilities::unpack<T>(temp, false);
+ }
+#endif
(void)comm;
+ (void)combiner;
return vec;
-#endif
+ }
+
+
+ template <typename T>
+ std::vector<T>
+ compute_set_union(const std::vector<T> &vec, const MPI_Comm &comm)
+ {
+ return Utilities::MPI::all_reduce<std::vector<T>>(
+ vec, comm, [](const auto &set_1, const auto &set_2) {
+ // merge vectors, sort, and eliminate duplicates -> mimic std::set<T>
+ auto result = set_1;
+ result.insert(result.end(), set_2.begin(), set_2.end());
+ std::sort(result.begin(), result.end());
+ result.erase(std::unique(result.begin(), result.end()), result.end());
+ return result;
+ });
}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Test Utilities::MPI::all_reduce().
+
+#include <deal.II/base/mpi.h>
+
+#include "../tests.h"
+
+
+template <typename T>
+void
+check(const std::function<std::vector<T>(const std::vector<T> &,
+ const std::vector<T> &)> &fu)
+{
+ const auto result = Utilities::MPI::all_reduce(
+ std::vector<T>{Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)},
+ MPI_COMM_WORLD,
+ fu);
+
+ for (const auto r : result)
+ deallog << r << " ";
+ deallog << std::endl;
+}
+
+
+int
+main(int argc, char *argv[])
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+ MPILogInitAll log;
+
+ // compute min
+ check<unsigned int>([](const auto &a, const auto &b) {
+ return std::vector<unsigned int>{std::min(a[0], b[0])};
+ });
+
+ // compute max
+ check<unsigned int>([](const auto &a, const auto &b) {
+ return std::vector<unsigned int>{std::max(a[0], b[0])};
+ });
+
+ // compute sum
+ check<unsigned int>([](const auto &a, const auto &b) {
+ return std::vector<unsigned int>{a[0] + b[0]};
+ });
+
+ // perform gather all
+ check<unsigned int>([](const auto &a, const auto &b) {
+ auto result = a;
+ result.insert(result.end(), b.begin(), b.end());
+ return result;
+ });
+}