From: Maximilian Bergbauer Date: Thu, 18 Mar 2021 14:45:15 +0000 (+0100) Subject: Add broadcast function. X-Git-Tag: v9.3.0-rc1~291^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=d061a0e40a2ce66aaa9e082c3f2105316e0cefb3;p=dealii.git Add broadcast function. --- diff --git a/doc/news/changes/minor/20210318MaximilianBergbauer b/doc/news/changes/minor/20210318MaximilianBergbauer new file mode 100644 index 0000000000..fb78ab18de --- /dev/null +++ b/doc/news/changes/minor/20210318MaximilianBergbauer @@ -0,0 +1,3 @@ +New: Adds the new generalized MPI_Bcast function Utilities::MPI::broadcast() for arbitrary data types T. +
+(Maximilian Bergbauer, 2021/03/18) diff --git a/include/deal.II/base/mpi.h b/include/deal.II/base/mpi.h index cb63553c3b..6a281b3798 100644 --- a/include/deal.II/base/mpi.h +++ b/include/deal.II/base/mpi.h @@ -1071,6 +1071,28 @@ namespace Utilities const T & object_to_send, const unsigned int root_process = 0); + /** + * Sends an object @p object_to_send from the process @p root_process + * to all other processes. + * + * A generalization of the classic `MPI_Bcast` function that accepts + * arbitrary data types `T`, as long as Utilities::pack() (which in turn + * uses `boost::serialize`, see in Utilities::pack() for details) accepts + * `T` as an argument. + * + * @param[in] comm MPI communicator. + * @param[in] object_to_send An object to send to all processes. + * @param[in] root_process The process that sends the object to all + * processes. By default the process with rank 0 is the root process. + * + * @return Every process receives the object sent by the @p root_process. + */ + template + T + broadcast(const MPI_Comm & comm, + 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 @@ -1463,6 +1485,37 @@ namespace Utilities # endif } + template + T + broadcast(const MPI_Comm & comm, + const T & object_to_send, + const unsigned int root_process) + { +# ifndef DEAL_II_WITH_MPI + (void)comm; + (void)root_process; + return object_to_send; +# else + const auto n_procs = dealii::Utilities::MPI::n_mpi_processes(comm); + AssertIndexRange(root_process, n_procs); + (void)n_procs; + + std::vector buffer = Utilities::pack(object_to_send, false); + unsigned int buffer_size = buffer.size(); + + // Exchanging the size of buffer + int ierr = MPI_Bcast(&buffer_size, 1, MPI_UNSIGNED, root_process, comm); + AssertThrowMPI(ierr); + + buffer.resize(buffer_size); + ierr = + MPI_Bcast(buffer.data(), buffer_size, MPI_CHAR, root_process, comm); + AssertThrowMPI(ierr); + + return Utilities::unpack(buffer, false); +# endif + } + # ifdef DEAL_II_WITH_MPI template diff --git a/include/deal.II/base/mpi.templates.h b/include/deal.II/base/mpi.templates.h index 1e2c6d960c..15a1045d1d 100644 --- a/include/deal.II/base/mpi.templates.h +++ b/include/deal.II/base/mpi.templates.h @@ -593,13 +593,7 @@ namespace Utilities } // 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 Utilities::unpack(temp, false); + return Utilities::MPI::broadcast(comm, result); } #endif (void)comm; diff --git a/tests/mpi/broadcast_01.cc b/tests/mpi/broadcast_01.cc new file mode 100644 index 0000000000..db92acf0de --- /dev/null +++ b/tests/mpi/broadcast_01.cc @@ -0,0 +1,83 @@ +// --------------------------------------------------------------------- +// +// 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::broadcast(). + +#include + +#include "../tests.h" + +struct TestObject +{ + std::vector int_vector; + std::string string; + + template + void + serialize(Archive &ar, const unsigned int /*version*/) + { + ar &int_vector; + ar &string; + } + + friend std::ostream & + operator<<(std::ostream &os, const TestObject &test_object); +}; + +std::ostream & +operator<<(std::ostream &os, const TestObject &test_object) +{ + for (const auto &i : test_object.int_vector) + os << i << " "; + + os << test_object.string; + + return os; +} + +void +check(const std::vector &int_vector, const std::string &string) +{ + const auto my_proc = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + TestObject test_object; + + if (my_proc == 0) + { + test_object.int_vector = int_vector; + test_object.string = string; + } + + const auto result = Utilities::MPI::broadcast(MPI_COMM_WORLD, test_object); + + deallog << result << " "; + deallog << std::endl; +} + + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll log; + + // broadcast test + check({1, 2, 3}, "test"); + + // broadcast test + check({-1, 0, 1}, "i love democracy"); +} diff --git a/tests/mpi/broadcast_01.mpirun=1.output b/tests/mpi/broadcast_01.mpirun=1.output new file mode 100644 index 0000000000..8407359e5a --- /dev/null +++ b/tests/mpi/broadcast_01.mpirun=1.output @@ -0,0 +1,3 @@ + +DEAL:0::1 2 3 test +DEAL:0::-1 0 1 i love democracy diff --git a/tests/mpi/broadcast_01.mpirun=5.output b/tests/mpi/broadcast_01.mpirun=5.output new file mode 100644 index 0000000000..2d4a59a717 --- /dev/null +++ b/tests/mpi/broadcast_01.mpirun=5.output @@ -0,0 +1,19 @@ + +DEAL:0::1 2 3 test +DEAL:0::-1 0 1 i love democracy + +DEAL:1::1 2 3 test +DEAL:1::-1 0 1 i love democracy + + +DEAL:2::1 2 3 test +DEAL:2::-1 0 1 i love democracy + + +DEAL:3::1 2 3 test +DEAL:3::-1 0 1 i love democracy + + +DEAL:4::1 2 3 test +DEAL:4::-1 0 1 i love democracy +