From: Peter Munch Date: Tue, 16 Feb 2021 09:34:30 +0000 (+0100) Subject: Introduce Utilities::MPI::allreduce X-Git-Tag: v9.3.0-rc1~404^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f6076952713ac7f9225120d84b28b96cdf67c4ca;p=dealii.git Introduce Utilities::MPI::allreduce --- diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index b138d3433f..d673e2cf8f 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -1008,6 +1008,21 @@ namespace Utilities const T & object_to_send, const unsigned int root_process = 0); + /** + * A function that combines values @p local_value from all processes + * via a user-specified binary operation @p combiner and distributes the + * result back to all processes. As such this function is similar to + * MPI_Allreduce (if it were implemented by a global reduction followed + * by a broadcast step) but due to the user-specified binary operation also + * general object types, including ones that store variable amounts of data, + * can be handled. + */ + template + T + all_reduce(const T & local_value, + const MPI_Comm & comm, + const std::function &combiner); + /** * Given a partitioned index set space, compute the owning MPI process rank * of each element of a second index set according to the partitioned index diff --git a/include/deal.II/base/mpi.templates.h b/include/deal.II/base/mpi.templates.h index e30dee3e3d..fa3cf3a31d 100644 --- a/include/deal.II/base/mpi.templates.h +++ b/include/deal.II/base/mpi.templates.h @@ -457,80 +457,98 @@ namespace Utilities template - std::vector - compute_set_union(const std::vector &vec, const MPI_Comm &comm) + T + all_reduce(const T & vec, + const MPI_Comm & comm, + const std::function &combiner) { #ifdef DEAL_II_WITH_MPI - // 1) collect vector entries and create union - std::vector 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 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(temp, false)); + } + else if (rank_send == rank) // process sends data + { + const std::vector 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 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(temp, false); + } +#endif (void)comm; + (void)combiner; return vec; -#endif + } + + + template + std::vector + compute_set_union(const std::vector &vec, const MPI_Comm &comm) + { + return Utilities::MPI::all_reduce>( + vec, comm, [](const auto &set_1, const auto &set_2) { + // merge vectors, sort, and eliminate duplicates -> mimic std::set + 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; + }); } diff --git a/source/base/mpi.inst.in b/source/base/mpi.inst.in index 9e1de666fe..0b59283041 100644 --- a/source/base/mpi.inst.in +++ b/source/base/mpi.inst.in @@ -67,6 +67,17 @@ for (S : MPI_SCALARS) const MPI_Comm &, const ArrayView &); + template S all_reduce( + const S & vec, + const MPI_Comm & comm, + const std::function &process); + + template std::vector all_reduce( + const std::vector & vec, + const MPI_Comm & comm, + const std::function(const std::vector &, + const std::vector &)> &process); + // The fixed-length array (i.e., things declared like T(&values)[N]) // versions of the functions above live in the header file mpi.h since the // length (N) is a compile-time constant. Those functions all call diff --git a/tests/mpi/allreduce_01.cc b/tests/mpi/allreduce_01.cc new file mode 100644 index 0000000000..46880fcef3 --- /dev/null +++ b/tests/mpi/allreduce_01.cc @@ -0,0 +1,68 @@ +// --------------------------------------------------------------------- +// +// 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 + +#include "../tests.h" + + +template +void +check(const std::function(const std::vector &, + const std::vector &)> &fu) +{ + const auto result = Utilities::MPI::all_reduce( + std::vector{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([](const auto &a, const auto &b) { + return std::vector{std::min(a[0], b[0])}; + }); + + // compute max + check([](const auto &a, const auto &b) { + return std::vector{std::max(a[0], b[0])}; + }); + + // compute sum + check([](const auto &a, const auto &b) { + return std::vector{a[0] + b[0]}; + }); + + // perform gather all + check([](const auto &a, const auto &b) { + auto result = a; + result.insert(result.end(), b.begin(), b.end()); + return result; + }); +} diff --git a/tests/mpi/allreduce_01.mpirun=1.output b/tests/mpi/allreduce_01.mpirun=1.output new file mode 100644 index 0000000000..8a9e3e0734 --- /dev/null +++ b/tests/mpi/allreduce_01.mpirun=1.output @@ -0,0 +1,5 @@ + +DEAL:0::0 +DEAL:0::0 +DEAL:0::0 +DEAL:0::0 diff --git a/tests/mpi/allreduce_01.mpirun=5.output b/tests/mpi/allreduce_01.mpirun=5.output new file mode 100644 index 0000000000..0b8dac6bdf --- /dev/null +++ b/tests/mpi/allreduce_01.mpirun=5.output @@ -0,0 +1,29 @@ + +DEAL:0::0 +DEAL:0::4 +DEAL:0::10 +DEAL:0::0 1 2 3 4 + +DEAL:1::0 +DEAL:1::4 +DEAL:1::10 +DEAL:1::0 1 2 3 4 + + +DEAL:2::0 +DEAL:2::4 +DEAL:2::10 +DEAL:2::0 1 2 3 4 + + +DEAL:3::0 +DEAL:3::4 +DEAL:3::10 +DEAL:3::0 1 2 3 4 + + +DEAL:4::0 +DEAL:4::4 +DEAL:4::10 +DEAL:4::0 1 2 3 4 +