]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add broadcast function. 11940/head
authorMaximilian Bergbauer <bergbauer@lnm.mw.tum.de>
Thu, 18 Mar 2021 14:45:15 +0000 (15:45 +0100)
committerMaximilian Bergbauer <bergbauer@lnm.mw.tum.de>
Mon, 22 Mar 2021 10:35:09 +0000 (11:35 +0100)
doc/news/changes/minor/20210318MaximilianBergbauer [new file with mode: 0644]
include/deal.II/base/mpi.h
include/deal.II/base/mpi.templates.h
tests/mpi/broadcast_01.cc [new file with mode: 0644]
tests/mpi/broadcast_01.mpirun=1.output [new file with mode: 0644]
tests/mpi/broadcast_01.mpirun=5.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20210318MaximilianBergbauer b/doc/news/changes/minor/20210318MaximilianBergbauer
new file mode 100644 (file)
index 0000000..fb78ab1
--- /dev/null
@@ -0,0 +1,3 @@
+New: Adds the new generalized MPI_Bcast function Utilities::MPI::broadcast() for arbitrary data types T.
+<br>
+(Maximilian Bergbauer, 2021/03/18)
index cb63553c3b339b7f5b905d7de33d4b851d3c94b3..6a281b3798fd63eddc4e82445d24613171be2e60 100644 (file)
@@ -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 <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
@@ -1463,6 +1485,37 @@ namespace Utilities
 #  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>
index 1e2c6d960cb687d62eb12e47e112c11df4a701ed..15a1045d1d4837e85a5ed46fa539ecfb5fd56d43 100644 (file)
@@ -593,13 +593,7 @@ namespace Utilities
             }
 
           // 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;
diff --git a/tests/mpi/broadcast_01.cc b/tests/mpi/broadcast_01.cc
new file mode 100644 (file)
index 0000000..db92acf
--- /dev/null
@@ -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 <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");
+}
diff --git a/tests/mpi/broadcast_01.mpirun=1.output b/tests/mpi/broadcast_01.mpirun=1.output
new file mode 100644 (file)
index 0000000..8407359
--- /dev/null
@@ -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 (file)
index 0000000..2d4a59a
--- /dev/null
@@ -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
+

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.