]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MPI_Allreduce is apparently unable to accept input and output arguments at the same...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 29 Sep 2014 15:12:57 +0000 (10:12 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 29 Sep 2014 15:12:57 +0000 (10:12 -0500)
Fix this by determining when &input=&output and in that case pass MPI_IN_PLACE.

14 files changed:
doc/news/changes.h
include/deal.II/base/mpi.h
tests/mpi/collective_02_array_in_place.cc [new file with mode: 0644]
tests/mpi/collective_02_array_in_place.mpirun=1.output [new file with mode: 0644]
tests/mpi/collective_02_array_in_place.mpirun=10.output [new file with mode: 0644]
tests/mpi/collective_02_array_in_place.mpirun=4.output [new file with mode: 0644]
tests/mpi/collective_02_vector.cc [new file with mode: 0644]
tests/mpi/collective_02_vector.mpirun=1.output [new file with mode: 0644]
tests/mpi/collective_02_vector.mpirun=10.output [new file with mode: 0644]
tests/mpi/collective_02_vector.mpirun=4.output [new file with mode: 0644]
tests/mpi/collective_02_vector_in_place.cc [new file with mode: 0644]
tests/mpi/collective_02_vector_in_place.mpirun=1.output [new file with mode: 0644]
tests/mpi/collective_02_vector_in_place.mpirun=10.output [new file with mode: 0644]
tests/mpi/collective_02_vector_in_place.mpirun=4.output [new file with mode: 0644]

index 4bb5c9fc1b75067f2d03a82d75389d9f93da2efd..62fbcd1efd4cbb980ab092c752a8108394bef710 100644 (file)
@@ -301,6 +301,13 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> Fixed: The vector and array versions of Utilities::MPI::sum() and
+  Utilities::MPI::max() produced segmentation faults with some MPI implementations
+  if the input and output arguments were the same. This is now fixed.
+  <br>
+  (Wolfgang Bangerth, 2014/09/29)
+  </li>
+
   <li> Fixed: Trying to have FE_Q(p) and FE_DGQ(r) elements next to each
   other in an hp::DoFHandler object led to assertions saying that these two
   elements don't know how to compute interface constraints where such
index 5bbac5053978829e482b3f88635a3f333e2e0d8c..a6943e36910a9f8c81e7ccfca55988cb70cf0d22 100644 (file)
@@ -138,6 +138,8 @@ namespace Utilities
      * array of length N. In other words, the i-th element of the results
      * array is the sum over the i-th entries of the input arrays from each
      * processor.
+     *
+     * Input and output arrays may be the same.
      */
     template <typename T, unsigned int N>
     inline
@@ -149,6 +151,8 @@ namespace Utilities
      * Like the previous function, but take the sums over the elements of a
      * std::vector. In other words, the i-th element of the results array is
      * the sum over the i-th entries of the input arrays from each processor.
+     *
+     * Input and output vectors may be the same.
      */
     template <typename T>
     inline
@@ -184,6 +188,8 @@ namespace Utilities
      * array of length N. In other words, the i-th element of the results
      * array is the maximum of the i-th entries of the input arrays from each
      * processor.
+     *
+     * Input and output arrays may be the same.
      */
     template <typename T, unsigned int N>
     inline
@@ -196,6 +202,8 @@ namespace Utilities
      * std::vector. In other words, the i-th element of the results array is
      * the maximum over the i-th entries of the input arrays from each
      * processor.
+     *
+     * Input and output vectors may be the same.
      */
     template <typename T>
     inline
@@ -415,7 +423,11 @@ namespace Utilities
               T (&sums)[N])
     {
 #ifdef DEAL_II_WITH_MPI
-      MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&values[0])),
+      MPI_Allreduce ((&values[0] != &sums[0]
+                      ?
+                      const_cast<void *>(static_cast<const void *>(&values[0]))
+                      :
+                      MPI_IN_PLACE),
                      &sums[0], N, internal::mpi_type_id(values), MPI_SUM,
                      mpi_communicator);
 #else
@@ -434,7 +446,11 @@ namespace Utilities
     {
 #ifdef DEAL_II_WITH_MPI
       sums.resize (values.size());
-      MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&values[0])),
+      MPI_Allreduce ((&values[0] != &sums[0]
+                      ?
+                      const_cast<void *>(static_cast<const void *>(&values[0]))
+                      :
+                      MPI_IN_PLACE),
                      &sums[0], values.size(), internal::mpi_type_id((T *)0), MPI_SUM,
                      mpi_communicator);
 #else
@@ -469,7 +485,11 @@ namespace Utilities
               T (&maxima)[N])
     {
 #ifdef DEAL_II_WITH_MPI
-      MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&values[0])),
+      MPI_Allreduce ((&values[0] != &maxima[0]
+                      ?
+                      const_cast<void *>(static_cast<const void *>(&values[0]))
+                      :
+                      MPI_IN_PLACE),
                      &maxima[0], N, internal::mpi_type_id(values), MPI_MAX,
                      mpi_communicator);
 #else
@@ -488,7 +508,11 @@ namespace Utilities
     {
 #ifdef DEAL_II_WITH_MPI
       maxima.resize (values.size());
-      MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&values[0])),
+      MPI_Allreduce ((&values[0] != &maxima[0]
+                      ?
+                      const_cast<void *>(static_cast<const void *>(&values[0]))
+                      :
+                      MPI_IN_PLACE),
                      &maxima[0], values.size(), internal::mpi_type_id((T *)0), MPI_MAX,
                      mpi_communicator);
 #else
diff --git a/tests/mpi/collective_02_array_in_place.cc b/tests/mpi/collective_02_array_in_place.cc
new file mode 100644 (file)
index 0000000..49f3005
--- /dev/null
@@ -0,0 +1,66 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2014 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check Utilities::MPI::sum() for arrays, but with input=output
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+#include <fstream>
+
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  unsigned int sums[2] = { 1, 2 };
+  Utilities::MPI::sum (sums,
+                       MPI_COMM_WORLD,
+                       sums);
+  Assert (sums[0] == numprocs, ExcInternalError());
+  Assert (sums[1] == 2*numprocs, ExcInternalError());
+
+  if (myid==0)
+    deallog << sums[0] << ' ' << sums[1] << std::endl;
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_WITH_MPI
+  Utilities::MPI::MPI_InitFinalize mpi (argc, argv);
+#else
+  (void)argc;
+  (void)argv;
+  compile_time_error;
+
+#endif
+
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("mpi");
+      test();
+      deallog.pop();
+    }
+  else
+    test();
+}
diff --git a/tests/mpi/collective_02_array_in_place.mpirun=1.output b/tests/mpi/collective_02_array_in_place.mpirun=1.output
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2
diff --git a/tests/mpi/collective_02_array_in_place.mpirun=10.output b/tests/mpi/collective_02_array_in_place.mpirun=10.output
new file mode 100644 (file)
index 0000000..97f11ae
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::10 20
diff --git a/tests/mpi/collective_02_array_in_place.mpirun=4.output b/tests/mpi/collective_02_array_in_place.mpirun=4.output
new file mode 100644 (file)
index 0000000..ad500b9
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::4 8
diff --git a/tests/mpi/collective_02_vector.cc b/tests/mpi/collective_02_vector.cc
new file mode 100644 (file)
index 0000000..f9125af
--- /dev/null
@@ -0,0 +1,68 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2014 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check Utilities::MPI::sum() for vectors
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+#include <fstream>
+
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  unsigned int values_[2] = { 1, 2 };
+  std::vector<unsigned int> values(&values_[0], &values_[2]);
+  std::vector<unsigned int> sums(2);
+  Utilities::MPI::sum (values,
+                       MPI_COMM_WORLD,
+                       sums);
+  Assert (sums[0] == numprocs, ExcInternalError());
+  Assert (sums[1] == 2*numprocs, ExcInternalError());
+
+  if (myid==0)
+    deallog << sums[0] << ' ' << sums[1] << std::endl;
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_WITH_MPI
+  Utilities::MPI::MPI_InitFinalize mpi (argc, argv);
+#else
+  (void)argc;
+  (void)argv;
+  compile_time_error;
+
+#endif
+
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("mpi");
+      test();
+      deallog.pop();
+    }
+  else
+    test();
+}
diff --git a/tests/mpi/collective_02_vector.mpirun=1.output b/tests/mpi/collective_02_vector.mpirun=1.output
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2
diff --git a/tests/mpi/collective_02_vector.mpirun=10.output b/tests/mpi/collective_02_vector.mpirun=10.output
new file mode 100644 (file)
index 0000000..97f11ae
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::10 20
diff --git a/tests/mpi/collective_02_vector.mpirun=4.output b/tests/mpi/collective_02_vector.mpirun=4.output
new file mode 100644 (file)
index 0000000..ad500b9
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::4 8
diff --git a/tests/mpi/collective_02_vector_in_place.cc b/tests/mpi/collective_02_vector_in_place.cc
new file mode 100644 (file)
index 0000000..d56c029
--- /dev/null
@@ -0,0 +1,67 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2014 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check Utilities::MPI::sum() for vectors, but with input=output
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/utilities.h>
+#include <fstream>
+
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  unsigned int values_[2] = { 1, 2 };
+  std::vector<unsigned int> sums(&values_[0], &values_[2]);
+  Utilities::MPI::sum (sums,
+                       MPI_COMM_WORLD,
+                       sums);
+  Assert (sums[0] == numprocs, ExcInternalError());
+  Assert (sums[1] == 2*numprocs, ExcInternalError());
+
+  if (myid==0)
+    deallog << sums[0] << ' ' << sums[1] << std::endl;
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_WITH_MPI
+  Utilities::MPI::MPI_InitFinalize mpi (argc, argv);
+#else
+  (void)argc;
+  (void)argv;
+  compile_time_error;
+
+#endif
+
+  if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("mpi");
+      test();
+      deallog.pop();
+    }
+  else
+    test();
+}
diff --git a/tests/mpi/collective_02_vector_in_place.mpirun=1.output b/tests/mpi/collective_02_vector_in_place.mpirun=1.output
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2
diff --git a/tests/mpi/collective_02_vector_in_place.mpirun=10.output b/tests/mpi/collective_02_vector_in_place.mpirun=10.output
new file mode 100644 (file)
index 0000000..97f11ae
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::10 20
diff --git a/tests/mpi/collective_02_vector_in_place.mpirun=4.output b/tests/mpi/collective_02_vector_in_place.mpirun=4.output
new file mode 100644 (file)
index 0000000..ad500b9
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::4 8

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.