]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add Utilities::MPI::compute_union
authorPeter Munch <peterrmuench@gmail.com>
Sat, 4 Apr 2020 18:18:24 +0000 (20:18 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 7 Apr 2020 08:23:05 +0000 (10:23 +0200)
doc/news/changes/minor/20200404PeterMunch [new file with mode: 0644]
include/deal.II/base/mpi.h
include/deal.II/base/mpi.templates.h
include/deal.II/base/mpi_tags.h
source/base/mpi.cc
tests/mpi/union_01.cc [new file with mode: 0644]
tests/mpi/union_01.mpirun=2.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20200404PeterMunch b/doc/news/changes/minor/20200404PeterMunch
new file mode 100644 (file)
index 0000000..20b83cd
--- /dev/null
@@ -0,0 +1,4 @@
+New: Add new function Utilities::MPI::compute_union(), which computes the union of given input 
+vectors of all processes in a given MPI communicator.
+<br>
+(Peter Munch, 2020/04/04)
index 6431e56ea15de8c1d590be5c03a74a264790d2ab..8c118a978bef75ac32a5df0baf178d67fc2b1309 100644 (file)
@@ -1024,6 +1024,17 @@ namespace Utilities
                         const IndexSet &indices_to_look_up,
                         const MPI_Comm &comm);
 
+    /**
+     * Compute the union of the input vectors @p vec of all processes in the
+     *   MPI communicator @p comm.
+     *
+     * @note This is a collective operation. The result will available on all
+     *   processes.
+     */
+    template <typename T>
+    std::vector<T>
+    compute_union(const std::vector<T> &vec, const MPI_Comm &comm);
+
 #ifndef DOXYGEN
     // declaration for an internal function that lives in mpi.templates.h
     namespace internal
index 4c0ea2c82121d43551e4b2da1fd735c5c91c415f..d99d6c484181e57d7d79a5a5fdaa36f7af0288a3 100644 (file)
@@ -28,6 +28,7 @@
 #include <deal.II/lac/sparse_matrix.h>
 #include <deal.II/lac/vector.h>
 
+#include <set>
 #include <vector>
 
 DEAL_II_NAMESPACE_OPEN
@@ -413,6 +414,76 @@ namespace Utilities
       internal::all_reduce(MPI_MIN, values, mpi_communicator, minima);
     }
 
+
+
+    template <typename T>
+    std::vector<T>
+    compute_union(const std::vector<T> &vec, const MPI_Comm &comm)
+    {
+#ifdef DEAL_II_WITH_MPI
+      // 1) collect vector entries and create union
+      std::vector<T> result = vec;
+
+      if (this_mpi_process(comm) == 0)
+        {
+          for (unsigned int i = 1; i < n_mpi_processes(comm); i++)
+            {
+              MPI_Status status;
+              const auto ierr_1 = MPI_Probe(MPI_ANY_SOURCE,
+                                            internal::Tags::compute_union,
+                                            comm,
+                                            &status);
+              AssertThrowMPI(ierr_1);
+
+              int        amount;
+              const auto ierr_2 =
+                MPI_Get_count(&status,
+                              internal::mpi_type_id(vec.data()),
+                              &amount);
+              AssertThrowMPI(ierr_2);
+
+              const unsigned int old_size = result.size();
+              result.resize(old_size + amount);
+
+              const auto ierr_3 = MPI_Recv(result.data() + old_size,
+                                           amount,
+                                           MPI_INT,
+                                           status.MPI_SOURCE,
+                                           status.MPI_TAG,
+                                           comm,
+                                           MPI_STATUS_IGNORE);
+              AssertThrowMPI(ierr_3);
+
+              std::sort(result.begin(), result.end());
+              result.erase(std::unique(result.begin(), result.end()),
+                           result.end());
+            }
+        }
+      else
+        {
+          const auto ierr = MPI_Send(vec.data(),
+                                     vec.size(),
+                                     internal::mpi_type_id(vec.data()),
+                                     0,
+                                     internal::Tags::compute_union,
+                                     comm);
+          AssertThrowMPI(ierr);
+        }
+
+      // 2) broadcast result
+      int size = result.size();
+      MPI_Bcast(&size, 1, MPI_INT, 0, comm);
+      result.resize(size);
+      MPI_Bcast(
+        result.data(), size, internal::mpi_type_id(vec.data()), 0, comm);
+
+      return result;
+#else
+      (void)comm;
+      return vec;
+#endif
+    }
+
   } // end of namespace MPI
 } // end of namespace Utilities
 
index 548df86b964efae707b7074b61a31ad09a8c62dc..4c5d2f533326c172fad3012aaa6e90fcfa2ca7a8 100644 (file)
@@ -132,6 +132,9 @@ namespace Utilities
           /// NoncontiguousPartitioner::update_values
           noncontiguous_partitioner_update_values,
 
+          // Utilities::MPI::compute_union
+          compute_union,
+
         };
       } // namespace Tags
     }   // namespace internal
index d0a62bc9bb3f0060148e57a92bb78ae6a58a3ab3..d6784fe361f4e91652dee2349041175ead6779b0 100644 (file)
@@ -1126,6 +1126,10 @@ namespace Utilities
       locked = false;
     }
 
+
+    template std::vector<unsigned int>
+    compute_union(const std::vector<unsigned int> &vec, const MPI_Comm &comm);
+
 #include "mpi.inst"
   } // end of namespace MPI
 } // end of namespace Utilities
diff --git a/tests/mpi/union_01.cc b/tests/mpi/union_01.cc
new file mode 100644 (file)
index 0000000..884e113
--- /dev/null
@@ -0,0 +1,51 @@
+// ---------------------------------------------------------------------
+//
+// 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Create a distributed triangulation with multigrid levels and copy it.
+
+#include <deal.II/base/mpi.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+void
+test(const MPI_Comm comm)
+{
+  std::vector<unsigned int> vector;
+
+  if (Utilities::MPI::this_mpi_process(comm) == 0)
+    vector = std::vector<unsigned int>{2, 5, 7};
+  else
+    vector = std::vector<unsigned int>{5, 8, 9, 10};
+
+  const auto result = Utilities::MPI::compute_union(vector, comm);
+
+  for (auto i : result)
+    deallog << i << " ";
+  deallog << std::endl;
+}
+
+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/mpi/union_01.mpirun=2.output b/tests/mpi/union_01.mpirun=2.output
new file mode 100644 (file)
index 0000000..ba6eac9
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL:0::2 5 7 8 9 10 
+
+DEAL:1::2 5 7 8 9 10 
+

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.