]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Add Utilities::MPI::max.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 26 Sep 2011 17:54:20 +0000 (17:54 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 26 Sep 2011 17:54:20 +0000 (17:54 +0000)
git-svn-id: https://svn.dealii.org/trunk@24433 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/base/utilities.h
tests/mpi/collective_02.cc
tests/mpi/collective_03.cc [new file with mode: 0644]
tests/mpi/collective_03/ncpu_1/cmp/generic [new file with mode: 0644]
tests/mpi/collective_03/ncpu_10/cmp/generic [new file with mode: 0644]
tests/mpi/collective_03/ncpu_4/cmp/generic [new file with mode: 0644]

index e030ab9e507110ebc69344c2ea5318352d8d44b1..5ea9804271ff4aa337b57d8af33837633f2d4b5d 100644 (file)
@@ -336,10 +336,8 @@ its "local" index or numbers::invalid_unsigned_int.
 <br>
 (Guido Kanschat, 2011/09/26)
 
-<li> New: There is now a function GridTools::volume() computing the volume
-
-<li> New: The function Utilities::MPI::sum
-function can be used to compute the sum of values over a number of
+<li> New: The functions Utilities::MPI::sum and Utilities::MPI::max
+function can be used to compute the sum or maximum of values over a number of
 MPI processes without having to deal with the underlying MPI functions.
 <br>
 (Wolfgang Bangerth, 2011/09/25)
@@ -366,7 +364,7 @@ threads kept pouncing on it. This is now fixed.
 <br>
 (Patrick Sodr&eacute;, 2011/09/07)
 
-<li> New: There is now a function GridTools::volume computing the volume
+<li> New: There is now a function GridTools::volume() computing the volume
 of a triangulation.
 <br>
 (Wolfgang Bangerth, 2011/09/02)
index 803d62e78fd808785aca2400f61176c05bc1c2e4..9b16ce443410be27b75d1713f1e9c9d331dd5206 100644 (file)
@@ -381,6 +381,22 @@ namespace Utilities
      */
     template <typename T>
     T sum (const T &t,
+          const MPI_Comm &mpi_communicator);
+
+    /**
+     * Return the maximum over all processors of the value @p t. This function
+     * is collective over all processors given in the communicator. If
+     * deal.II is not configured for use of MPI, this function simply
+     * returns the value of @p t. This function corresponds to the
+     * <code>MPI_Allreduce</code> function, i.e. all processors receive
+     * the result of this operation.
+     *
+     * @note This function is only implemented for certain template
+     * arguments <code>T</code>, namely <code>float, double, int,
+     * unsigned int</code>.
+     */
+    template <typename T>
+    T max (const T &t,
           const MPI_Comm &mpi_communicator);
 
                                     /**
@@ -983,6 +999,7 @@ namespace Utilities
 #endif
     }
 
+
     template <typename T>
     inline
     T sum (const T &t,
@@ -999,6 +1016,24 @@ namespace Utilities
       return t;
 #endif
     }
+
+
+    template <typename T>
+    inline
+    T max (const T &t,
+          const MPI_Comm &mpi_communicator)
+    {
+#ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+      T sum;
+      MPI_Allreduce (const_cast<void*>(static_cast<const void*>(&t)),
+                    &sum, 1, internal::mpi_type_id(&t), MPI_MAX,
+                    mpi_communicator);
+      return sum;
+#else
+      (void)mpi_communicator;
+      return t;
+#endif
+    }
   }
 }
 
index d053862530c3e4b73cdb8ddabdba065c00f36837..4f1cc1ea932d23775cb34f58d81efd5309c86d6a 100644 (file)
@@ -12,7 +12,7 @@
 //---------------------------------------------------------------------------
 
 
-// check Utilities::System::calculate_collective_mpi_sum()
+// check Utilities::MPI::sum()
 
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
diff --git a/tests/mpi/collective_03.cc b/tests/mpi/collective_03.cc
new file mode 100644 (file)
index 0000000..89a6bee
--- /dev/null
@@ -0,0 +1,74 @@
+//---------------------------------------------------------------------------
+//    $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()
+
+#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);
+
+  int int_sum, uint_sum, double_sum, float_sum;
+  int_sum
+    =
+    Utilities::MPI::max<int>(myid+1,
+                            MPI_COMM_WORLD);
+  uint_sum
+    =
+    Utilities::MPI::max<unsigned int>(myid+1,
+                                     MPI_COMM_WORLD);
+  float_sum
+    =
+    Utilities::MPI::max<float>(myid+1,
+                              MPI_COMM_WORLD);
+  double_sum
+    =
+    Utilities::MPI::max<double>(myid+1,
+                               MPI_COMM_WORLD);
+
+  if (myid==0)
+    deallog << int_sum << ' ' << uint_sum << ' ' << double_sum << ' ' << float_sum << 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").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/ncpu_1/cmp/generic b/tests/mpi/collective_03/ncpu_1/cmp/generic
new file mode 100644 (file)
index 0000000..4723d36
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::1 1 1 1
diff --git a/tests/mpi/collective_03/ncpu_10/cmp/generic b/tests/mpi/collective_03/ncpu_10/cmp/generic
new file mode 100644 (file)
index 0000000..c72ae8e
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::10 10 10 10
diff --git a/tests/mpi/collective_03/ncpu_4/cmp/generic b/tests/mpi/collective_03/ncpu_4/cmp/generic
new file mode 100644 (file)
index 0000000..b0dbea0
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL:mpi::4 4 4 4

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.