]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce Utilities::MPI::mpi_processes_within_communicator 11137/head
authorPeter Munch <peterrmuench@gmail.com>
Tue, 3 Nov 2020 13:44:23 +0000 (14:44 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 3 Nov 2020 13:44:23 +0000 (14:44 +0100)
include/deal.II/base/mpi.h
include/deal.II/matrix_free/matrix_free.templates.h
source/base/mpi.cc

index 19ae43725f47e871b95c7d0d230451983d5af57c..f04914626c1c8b984fa69b258688b007d3eb54b5 100644 (file)
@@ -152,6 +152,14 @@ namespace Utilities
     unsigned int
     this_mpi_process(const MPI_Comm &mpi_communicator);
 
+    /**
+     * Return a vector of the ranks (within @p comm_large) of a subset of
+     * processes specified by @p comm_small.
+     */
+    const std::vector<unsigned int>
+    mpi_processes_within_communicator(const MPI_Comm &comm_large,
+                                      const MPI_Comm &comm_small);
+
     /**
      * Consider an unstructured communication pattern where every process in
      * an MPI universe wants to send some data to a subset of the other
index 6933925da7fecf23f8e2c00fb78f4e6ad785de51..f7b58e75266bbaeb4aa9322627fd375ecb206b35 100644 (file)
@@ -1495,26 +1495,9 @@ MatrixFree<dim, Number, VectorizedArrayType>::initialize_indices(
 
     if (is_non_buffering_sm_supported)
       {
-        // TODO: move into Utilities::MPI
-        const auto procs_of_sm = [](const MPI_Comm &comm,
-                                    const MPI_Comm &comm_shared) {
-          if (Utilities::MPI::job_supports_mpi() == false)
-            return std::vector<unsigned int>{0};
-
-          const unsigned int rank = Utilities::MPI::this_mpi_process(comm);
-          const unsigned int size =
-            Utilities::MPI::n_mpi_processes(comm_shared);
-
-          std::vector<unsigned int> ranks(size);
-          MPI_Allgather(
-            &rank, 1, MPI_UNSIGNED, ranks.data(), 1, MPI_UNSIGNED, comm_shared);
-
-          return ranks;
-        };
-
         // gather the ranks of the shared-memory domain
-        const std::vector<unsigned int> sm_procs =
-          procs_of_sm(task_info.communicator, communicator_sm);
+        const auto sm_procs = Utilities::MPI::mpi_processes_within_communicator(
+          task_info.communicator, communicator_sm);
 
         // get my rank within the shared-memory domain
         const unsigned int my_sm_pid = std::distance(
index 269aea2e59ea802e989de8f3e867bade734209ee..41b0a37c39ca45ad6de1d3d3cb98dd2567ef2b03 100644 (file)
@@ -135,6 +135,26 @@ namespace Utilities
     }
 
 
+
+    const std::vector<unsigned int>
+    mpi_processes_within_communicator(const MPI_Comm &comm_large,
+                                      const MPI_Comm &comm_small)
+    {
+      if (Utilities::MPI::job_supports_mpi() == false)
+        return std::vector<unsigned int>{0};
+
+      const unsigned int rank = Utilities::MPI::this_mpi_process(comm_large);
+      const unsigned int size = Utilities::MPI::n_mpi_processes(comm_small);
+
+      std::vector<unsigned int> ranks(size);
+      MPI_Allgather(
+        &rank, 1, MPI_UNSIGNED, ranks.data(), 1, MPI_UNSIGNED, comm_small);
+
+      return ranks;
+    }
+
+
+
     MPI_Comm
     duplicate_communicator(const MPI_Comm &mpi_communicator)
     {
@@ -686,6 +706,14 @@ namespace Utilities
 
 
 
+    const std::vector<unsigned int>
+    mpi_processes_within_communicator(const MPI_Comm &, const MPI_Comm &)
+    {
+      return std::vector<unsigned int>{0};
+    }
+
+
+
     std::vector<IndexSet>
     create_ascending_partitioning(const MPI_Comm & /*comm*/,
                                   const IndexSet::size_type local_size)

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.