]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add array versions of the sum/max operators.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 26 Sep 2011 18:26:56 +0000 (18:26 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 26 Sep 2011 18:26:56 +0000 (18:26 +0000)
git-svn-id: https://svn.dealii.org/trunk@24434 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/base/utilities.h
tests/mpi/collective_02_array.cc [new file with mode: 0644]
tests/mpi/collective_02_array/ncpu_1/cmp/generic [new file with mode: 0644]
tests/mpi/collective_02_array/ncpu_10/cmp/generic [new file with mode: 0644]
tests/mpi/collective_02_array/ncpu_4/cmp/generic [new file with mode: 0644]
tests/mpi/collective_03_array.cc [new file with mode: 0644]
tests/mpi/collective_03_array/ncpu_1/cmp/generic [new file with mode: 0644]
tests/mpi/collective_03_array/ncpu_10/cmp/generic [new file with mode: 0644]
tests/mpi/collective_03_array/ncpu_4/cmp/generic [new file with mode: 0644]

index 9b16ce443410be27b75d1713f1e9c9d331dd5206..f6e48b48fae15a0a0f9cf576260b8bb14814f083 100644 (file)
@@ -383,6 +383,22 @@ namespace Utilities
     T sum (const T &t,
           const MPI_Comm &mpi_communicator);
 
+                                    /**
+                                     * Like the previous function,
+                                     * but take the sums over the
+                                     * elements of an 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.
+                                     */
+    template <typename T, unsigned int N>
+    inline
+    void sum (const T (&values)[N],
+             const MPI_Comm &mpi_communicator,
+             T (&sums)[N]);
+
     /**
      * Return the maximum over all processors of the value @p t. This function
      * is collective over all processors given in the communicator. If
@@ -399,6 +415,22 @@ namespace Utilities
     T max (const T &t,
           const MPI_Comm &mpi_communicator);
 
+                                    /**
+                                     * Like the previous function,
+                                     * but take the maxima over the
+                                     * elements of an 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.
+                                     */
+    template <typename T, unsigned int N>
+    inline
+    void max (const T (&values)[N],
+             const MPI_Comm &mpi_communicator,
+             T (&maxima)[N]);
+
                                     /**
                                      * Data structure to store the result of
                                      * min_max_avg().
@@ -1018,6 +1050,24 @@ namespace Utilities
     }
 
 
+    template <typename T, unsigned int N>
+    inline
+    void sum (const T (&values)[N],
+             const MPI_Comm &mpi_communicator,
+             T (&sums)[N])
+    {
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+      MPI_Allreduce (const_cast<void*>(static_cast<const void*>(&values[0])),
+                    &sums[0], N, internal::mpi_type_id(values), MPI_SUM,
+                    mpi_communicator);
+#else
+      (void)mpi_communicator;
+      for (unsigned int i=0; i<N; ++i)
+       sums[i] = values[i];
+#endif
+    }
+
+
     template <typename T>
     inline
     T max (const T &t,
@@ -1034,6 +1084,24 @@ namespace Utilities
       return t;
 #endif
     }
+
+
+    template <typename T, unsigned int N>
+    inline
+    void max (const T (&values)[N],
+             const MPI_Comm &mpi_communicator,
+             T (&maxima)[N])
+    {
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+      MPI_Allreduce (const_cast<void*>(static_cast<const void*>(&values[0])),
+                    &maxima[0], N, internal::mpi_type_id(values), MPI_MAX,
+                    mpi_communicator);
+#else
+      (void)mpi_communicator;
+      for (unsigned int i=0; i<N; ++i)
+       maxima[i] = values[i];
+#endif
+    }
   }
 }
 
diff --git a/tests/mpi/collective_02_array.cc b/tests/mpi/collective_02_array.cc
new file mode 100644 (file)
index 0000000..35f9971
--- /dev/null
@@ -0,0 +1,64 @@
+//---------------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2009, 2010, 2011 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//---------------------------------------------------------------------------
+
+
+// check Utilities::MPI::sum() for arrays
+
+#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 };
+  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_COMPILER_SUPPORTS_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_file_for_mpi("collective_02_array").c_str());
+      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/ncpu_1/cmp/generic b/tests/mpi/collective_02_array/ncpu_1/cmp/generic
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/ncpu_10/cmp/generic b/tests/mpi/collective_02_array/ncpu_10/cmp/generic
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/ncpu_4/cmp/generic b/tests/mpi/collective_02_array/ncpu_4/cmp/generic
new file mode 100644 (file)
index 0000000..ad500b9
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::4 8
diff --git a/tests/mpi/collective_03_array.cc b/tests/mpi/collective_03_array.cc
new file mode 100644 (file)
index 0000000..b415514
--- /dev/null
@@ -0,0 +1,64 @@
+//---------------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2009, 2010, 2011 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//---------------------------------------------------------------------------
+
+
+// check Utilities::MPI::max() for arrays
+
+#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 };
+  unsigned int maxima[2];
+  Utilities::MPI::max (values,
+                      MPI_COMM_WORLD,
+                      maxima);
+  Assert (maxima[0] == 1, ExcInternalError());
+  Assert (maxima[1] == 2, ExcInternalError());
+
+  if (myid==0)
+    deallog << maxima[0] << ' ' << maxima[1] << std::endl;
+}
+
+
+int main(int argc, char *argv[])
+{
+#ifdef DEAL_II_COMPILER_SUPPORTS_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_file_for_mpi("collective_03_array").c_str());
+      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_03_array/ncpu_1/cmp/generic b/tests/mpi/collective_03_array/ncpu_1/cmp/generic
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2
diff --git a/tests/mpi/collective_03_array/ncpu_10/cmp/generic b/tests/mpi/collective_03_array/ncpu_10/cmp/generic
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2
diff --git a/tests/mpi/collective_03_array/ncpu_4/cmp/generic b/tests/mpi/collective_03_array/ncpu_4/cmp/generic
new file mode 100644 (file)
index 0000000..6e10373
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 2

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.