]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add AnonymousConsensusAlgorithmProcess 9574/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 28 Feb 2020 13:21:13 +0000 (14:21 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 2 Mar 2020 07:52:13 +0000 (08:52 +0100)
include/deal.II/base/mpi.h
include/deal.II/base/mpi_consensus_algorithm.h [new file with mode: 0644]
tests/base/consensus_algorithm_01.cc [new file with mode: 0644]
tests/base/consensus_algorithm_01.mpirun=4.output [new file with mode: 0644]

index 7e31c2f14e17beaf8abee99307796c05204914b3..ecb6ff11da11b30912864ea9c6bf9b6c3f8328c9 100644 (file)
@@ -1069,6 +1069,8 @@ namespace Utilities
                   const std::vector<T2> &recv_buffer);
     };
 
+
+
     /**
      * A base class for algorithms that implement the task of coming up with
      * communication patterns to retrieve data from other processes in a
diff --git a/include/deal.II/base/mpi_consensus_algorithm.h b/include/deal.II/base/mpi_consensus_algorithm.h
new file mode 100644 (file)
index 0000000..e9c13eb
--- /dev/null
@@ -0,0 +1,202 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_mpi_consensus_algorithm_h
+#define dealii_mpi_consensus_algorithm_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/mpi.templates.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+
+namespace Utilities
+{
+  namespace MPI
+  {
+    /**
+     * This class implements Utilities::MPI::ConsensusAlgorithmProcess,
+     * using user-provided function wrappers.
+     * The advantage of this class is that users do not have to write their
+     * own implementation but can register lambda functions directly.
+     */
+    template <typename T1, typename T2>
+    class AnonymousConsensusAlgorithmProcess
+      : public ConsensusAlgorithmProcess<T1, T2>
+    {
+    public:
+      /**
+       * Register functions that should be called for implementing the interface
+       * of ConsensusAlgorithmProcess.
+       *
+       * @param function_compute_targets called during `compute_targets`.
+       * @param function_create_request called during `create_request`.
+       * @param function_answer_request called during `answer_request`.
+       * @param function_prepare_buffer_for_answer called during
+       *   `prepare_buffer_for_answer`.
+       * @param function_read_answer called during `read_answer`.
+       */
+      AnonymousConsensusAlgorithmProcess(
+        const std::function<std::vector<unsigned int>()>
+          &function_compute_targets,
+        const std::function<void(const unsigned int, std::vector<T1> &)>
+          &function_create_request =
+            [](const unsigned int, std::vector<T1> &) {},
+        const std::function<void(const unsigned int,
+                                 const std::vector<T1> &,
+                                 std::vector<T2> &)> &function_answer_request =
+          [](const unsigned int, const std::vector<T1> &, std::vector<T2> &) {},
+        const std::function<void(const unsigned int, std::vector<T2> &)>
+          &function_prepare_buffer_for_answer =
+            [](const unsigned int, std::vector<T2> &) {},
+        const std::function<void(const unsigned int, const std::vector<T2> &)>
+          &function_read_answer =
+            [](const unsigned int, const std::vector<T2> &) {});
+
+      /**
+       * @copydoc ConsensusAlgorithmProcess::compute_targets()
+       */
+      std::vector<unsigned int>
+      compute_targets() override;
+
+      /**
+       * @copydoc ConsensusAlgorithmProcess::create_request()
+       */
+      void
+      create_request(const unsigned int other_rank,
+                     std::vector<T1> &  send_buffer) override;
+
+      /**
+       * @copydoc ConsensusAlgorithmProcess::answer_request()
+       */
+      void
+      answer_request(const unsigned int     other_rank,
+                     const std::vector<T1> &buffer_recv,
+                     std::vector<T2> &      request_buffer) override;
+
+      /**
+       * @copydoc ConsensusAlgorithmProcess::prepare_buffer_for_answer()
+       */
+      void
+      prepare_buffer_for_answer(const unsigned int other_rank,
+                                std::vector<T2> &  recv_buffer) override;
+
+      /**
+       * @copydoc ConsensusAlgorithmProcess::read_answer()
+       */
+      void
+      read_answer(const unsigned int     other_rank,
+                  const std::vector<T2> &recv_buffer) override;
+
+    private:
+      const std::function<std::vector<unsigned int>()> function_compute_targets;
+      const std::function<void(const int, std::vector<T1> &)>
+        function_create_request;
+      const std::function<
+        void(const unsigned int, const std::vector<T1> &, std::vector<T2> &)>
+        function_answer_request;
+      const std::function<void(const int, std::vector<T2> &)>
+        function_prepare_buffer_for_answer;
+      const std::function<void(const int, const std::vector<T2> &)>
+        function_read_answer;
+    };
+
+
+
+    template <typename T1, typename T2>
+    AnonymousConsensusAlgorithmProcess<T1, T2>::
+      AnonymousConsensusAlgorithmProcess(
+        const std::function<std::vector<unsigned int>()>
+          &function_compute_targets,
+        const std::function<void(const unsigned int, std::vector<T1> &)>
+          &                                           function_create_request,
+        const std::function<void(const unsigned int,
+                                 const std::vector<T1> &,
+                                 std::vector<T2> &)> &function_answer_request,
+        const std::function<void(const unsigned int, std::vector<T2> &)>
+          &function_prepare_buffer_for_answer,
+        const std::function<void(const unsigned int, const std::vector<T2> &)>
+          &function_read_answer)
+      : function_compute_targets(function_compute_targets)
+      , function_create_request(function_create_request)
+      , function_answer_request(function_answer_request)
+      , function_prepare_buffer_for_answer(function_prepare_buffer_for_answer)
+      , function_read_answer(function_read_answer)
+    {}
+
+
+
+    template <typename T1, typename T2>
+    std::vector<unsigned int>
+    AnonymousConsensusAlgorithmProcess<T1, T2>::compute_targets()
+    {
+      return function_compute_targets();
+    }
+
+
+
+    template <typename T1, typename T2>
+    void
+    AnonymousConsensusAlgorithmProcess<T1, T2>::create_request(
+      const unsigned int other_rank,
+      std::vector<T1> &  send_buffer)
+    {
+      function_create_request(other_rank, send_buffer);
+    }
+
+
+
+    template <typename T1, typename T2>
+    void
+    AnonymousConsensusAlgorithmProcess<T1, T2>::answer_request(
+      const unsigned int     other_rank,
+      const std::vector<T1> &buffer_recv,
+      std::vector<T2> &      request_buffer)
+    {
+      function_answer_request(other_rank, buffer_recv, request_buffer);
+    }
+
+
+
+    template <typename T1, typename T2>
+    void
+    AnonymousConsensusAlgorithmProcess<T1, T2>::prepare_buffer_for_answer(
+      const unsigned int other_rank,
+      std::vector<T2> &  recv_buffer)
+    {
+      function_prepare_buffer_for_answer(other_rank, recv_buffer);
+    }
+
+
+
+    template <typename T1, typename T2>
+    void
+    AnonymousConsensusAlgorithmProcess<T1, T2>::read_answer(
+      const unsigned int     other_rank,
+      const std::vector<T2> &recv_buffer)
+    {
+      function_read_answer(other_rank, recv_buffer);
+    }
+
+
+  } // end of namespace MPI
+} // end of namespace Utilities
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
diff --git a/tests/base/consensus_algorithm_01.cc b/tests/base/consensus_algorithm_01.cc
new file mode 100644 (file)
index 0000000..562bbe7
--- /dev/null
@@ -0,0 +1,71 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 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 AnonymousConsensusAlgorithmProcess.
+
+#include <deal.II/base/mpi_consensus_algorithm.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+void
+test(const MPI_Comm &comm)
+{
+  const unsigned int my_rank = dealii::Utilities::MPI::this_mpi_process(comm);
+  const unsigned int n_rank  = dealii::Utilities::MPI::n_mpi_processes(comm);
+
+  using T1 = unsigned int;
+  using T2 = unsigned int;
+
+  dealii::Utilities::MPI::AnonymousConsensusAlgorithmProcess<T1, T2> process(
+    [&]() {
+      std::vector<unsigned int> result{(my_rank + 1) % n_rank};
+      return result;
+    },
+    [&](const unsigned int other_rank, std::vector<T1> &send_buffer) {
+      send_buffer.push_back(my_rank);
+    },
+    [&](const unsigned int &   other_rank,
+        const std::vector<T1> &buffer_recv,
+        std::vector<T2> &      request_buffer) {
+      AssertDimension(other_rank, buffer_recv.front());
+      deallog << "ConsensusAlgorithmProcess::answer_request() passed!"
+              << std::endl;
+      request_buffer.push_back(my_rank);
+    },
+    [&](const unsigned int other_rank, std::vector<T2> &recv_buffer) {
+      recv_buffer.resize(1);
+    },
+    [&](const unsigned int other_rank, const std::vector<T2> &recv_buffer) {
+      AssertDimension(other_rank, recv_buffer.front());
+      deallog << "ConsensusAlgorithmProcess::function_read_answer() passed!"
+              << std::endl;
+    });
+  dealii::Utilities::MPI::ConsensusAlgorithmSelector<T1, T2>(process, comm)
+    .run();
+}
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  const MPI_Comm comm = MPI_COMM_WORLD;
+
+  test(comm);
+}
diff --git a/tests/base/consensus_algorithm_01.mpirun=4.output b/tests/base/consensus_algorithm_01.mpirun=4.output
new file mode 100644 (file)
index 0000000..8faf9d9
--- /dev/null
@@ -0,0 +1,15 @@
+
+DEAL:0::ConsensusAlgorithmProcess::answer_request() passed!
+DEAL:0::ConsensusAlgorithmProcess::function_read_answer() passed!
+
+DEAL:1::ConsensusAlgorithmProcess::answer_request() passed!
+DEAL:1::ConsensusAlgorithmProcess::function_read_answer() passed!
+
+
+DEAL:2::ConsensusAlgorithmProcess::answer_request() passed!
+DEAL:2::ConsensusAlgorithmProcess::function_read_answer() passed!
+
+
+DEAL:3::ConsensusAlgorithmProcess::answer_request() passed!
+DEAL:3::ConsensusAlgorithmProcess::function_read_answer() passed!
+

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.