--- /dev/null
+New: Adds the new generalized MPI_Bcast function Utilities::MPI::broadcast() for arbitrary data types T.
+<br>
+(Maximilian Bergbauer, 2021/03/18)
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 <typename T>
+ 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
# endif
}
+ template <typename T>
+ 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<char> 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<T>(buffer, false);
+# endif
+ }
+
# ifdef DEAL_II_WITH_MPI
template <class Iterator, typename Number>
}
// 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 Utilities::unpack<T>(temp, false);
+ return Utilities::MPI::broadcast(comm, result);
}
#endif
(void)comm;
--- /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::broadcast().
+
+#include <deal.II/base/mpi.h>
+
+#include "../tests.h"
+
+struct TestObject
+{
+ std::vector<int> int_vector;
+ std::string string;
+
+ template <class Archive>
+ 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> &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");
+}
--- /dev/null
+
+DEAL:0::1 2 3 test
+DEAL:0::-1 0 1 i love democracy
--- /dev/null
+
+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
+