]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move templates in mpi.h to mpi.templates.h.
authorDavid Wells <wellsd2@rpi.edu>
Mon, 30 May 2016 20:11:05 +0000 (16:11 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Mon, 30 May 2016 23:35:58 +0000 (19:35 -0400)
cmake/config/template-arguments.in
include/deal.II/base/mpi.h
include/deal.II/base/mpi.templates.h [new file with mode: 0644]
source/base/CMakeLists.txt
source/base/mpi.cc
source/base/mpi.inst.in [new file with mode: 0644]
tests/mpi/renumber_cuthill_mckee.cc
tests/mpi/renumber_cuthill_mckee_02.cc
tests/tests.h

index b3e4595a62c955a302d31e5277cc82fc8114e275..adf41b2ab8217e0860e1cb2c8c84f8f859f27d79 100644 (file)
@@ -2,6 +2,14 @@ BOOL            := { true; false }
 REAL_SCALARS    := { double; float }
 COMPLEX_SCALARS := { std::complex<double>;
                      std::complex<float> }
+MPI_SCALARS     := { int;
+                     long int;
+                     unsigned int;
+                     unsigned long int;
+                     unsigned long long int;
+                     float;
+                     double;
+                     long double }
 
 DERIVATIVE_TENSORS := { double;
                         Tensor<1,deal_II_dimension>;
index 943cb75f1fe426cac54abbf2135480b068ed4ec7..735c4b073105c6bc486da0e8a488ef10da54c6cc 100644 (file)
@@ -434,282 +434,43 @@ namespace Utilities
      */
     bool job_supports_mpi ();
 
+#ifndef DOXYGEN
+    // declaration for an internal function that lives in mpi.templates.h
     namespace internal
     {
-#ifdef DEAL_II_WITH_MPI
-      /**
-       * Return the corresponding MPI data type id for the argument given.
-       */
-      inline MPI_Datatype mpi_type_id (const int *)
-      {
-        return MPI_INT;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const long int *)
-      {
-        return MPI_LONG;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const unsigned int *)
-      {
-        return MPI_UNSIGNED;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const unsigned long int *)
-      {
-        return MPI_UNSIGNED_LONG;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const unsigned long long int *)
-      {
-        return MPI_UNSIGNED_LONG_LONG;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const float *)
-      {
-        return MPI_FLOAT;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const double *)
-      {
-        return MPI_DOUBLE;
-      }
-
-
-      inline MPI_Datatype mpi_type_id (const long double *)
-      {
-        return MPI_LONG_DOUBLE;
-      }
-#endif
-
       template <typename T>
-      inline
       void all_reduce (const MPI_Op      &mpi_op,
                        const T *const    values,
                        const MPI_Comm    &mpi_communicator,
                        T                 *output,
-                       const std::size_t  size)
-      {
-#ifdef DEAL_II_WITH_MPI
-        if (job_supports_mpi())
-          {
-            MPI_Allreduce (values != output
-                           ?
-                           // TODO This const_cast is only needed for older
-                           // (e.g., openMPI 1.6, released in 2012)
-                           // implementations of MPI-2. It is not needed as of
-                           // MPI-3 and we should remove it at some point in
-                           // the future.
-                           const_cast<void *>(static_cast<const void *> (values))
-                           :
-                           MPI_IN_PLACE,
-                           static_cast<void *>(output),
-                           static_cast<int>(size),
-                           internal::mpi_type_id(values),
-                           mpi_op,
-                           mpi_communicator);
-          }
-        else
-#endif
-          {
-            (void)mpi_op;
-            (void)mpi_communicator;
-            for (std::size_t i=0; i<size; ++i)
-              output[i] = values[i];
-          }
-      }
-
-      template <typename T>
-      inline
-      T all_reduce (const MPI_Op &mpi_op,
-                    const T &t,
-                    const MPI_Comm &mpi_communicator)
-      {
-        T output;
-        all_reduce(mpi_op, &t, mpi_communicator, &output, 1);
-        return output;
-      }
-
-      template <typename T, unsigned int N>
-      inline
-      void all_reduce (const MPI_Op &mpi_op,
-                       const T (&values)[N],
-                       const MPI_Comm &mpi_communicator,
-                       T (&output)[N])
-      {
-        all_reduce(mpi_op, values, mpi_communicator, output, N);
-      }
-
-      template <typename T>
-      inline
-      void all_reduce (const MPI_Op &mpi_op,
-                       const std::vector<T> &values,
-                       const MPI_Comm       &mpi_communicator,
-                       std::vector<T>       &output)
-      {
-        Assert(values.size() == output.size(),
-               ExcDimensionMismatch(values.size(), output.size()));
-        all_reduce(mpi_op, &values[0], mpi_communicator, &output[0], values.size());
-      }
-
-      template <typename T>
-      inline
-      void all_reduce (const MPI_Op    &mpi_op,
-                       const Vector<T> &values,
-                       const MPI_Comm  &mpi_communicator,
-                       Vector<T>  &output)
-      {
-        Assert(values.size() == output.size(),
-               ExcDimensionMismatch(values.size(), output.size()));
-        all_reduce(mpi_op, values.begin(), mpi_communicator, output.begin(), values.size());
-      }
+                       const std::size_t  size);
     }
 
-
-    template <typename T>
-    inline
-    T sum (const T &t,
-           const MPI_Comm &mpi_communicator)
-    {
-      return internal::all_reduce(MPI_SUM, t, mpi_communicator);
-    }
-
-
+    // Since these depend on N they must live in the header file
     template <typename T, unsigned int N>
-    inline
     void sum (const T (&values)[N],
               const MPI_Comm &mpi_communicator,
               T (&sums)[N])
     {
-      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums);
-    }
-
-
-    template <typename T>
-    inline
-    void sum (const std::vector<T> &values,
-              const MPI_Comm       &mpi_communicator,
-              std::vector<T>       &sums)
-    {
-      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums);
+      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums, N);
     }
 
-    template <typename T>
-    inline
-    void sum (const Vector<T> &values,
-              const MPI_Comm &mpi_communicator,
-              Vector<T> &sums)
-    {
-      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums);
-    }
-
-
-    template <int rank, int dim, typename Number>
-    inline
-    Tensor<rank,dim,Number>
-    sum (const Tensor<rank,dim,Number> &local,
-         const MPI_Comm &mpi_communicator)
-    {
-      const unsigned int n_entries = Tensor<rank,dim,Number>::n_independent_components;
-      Number entries[ Tensor<rank,dim,Number>::n_independent_components ];
-
-      for (unsigned int i=0; i< n_entries; ++i)
-        entries[i] = local[ local.unrolled_to_component_indices(i) ];
-
-      Number global_entries[ Tensor<rank,dim,Number>::n_independent_components ];
-      dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries );
-
-      Tensor<rank,dim,Number> global;
-      for (unsigned int i=0; i< n_entries; ++i)
-        global[ global.unrolled_to_component_indices(i) ] = global_entries[i];
-
-      return global;
-    }
-
-    template <int rank, int dim, typename Number>
-    inline
-    SymmetricTensor<rank,dim,Number>
-    sum (const SymmetricTensor<rank,dim,Number> &local,
-         const MPI_Comm &mpi_communicator)
-    {
-      const unsigned int n_entries = SymmetricTensor<rank,dim,Number>::n_independent_components;
-      Number entries[ SymmetricTensor<rank,dim,Number>::n_independent_components ];
-
-      for (unsigned int i=0; i< n_entries; ++i)
-        entries[i] = local[ local.unrolled_to_component_indices(i) ];
-
-      Number global_entries[ SymmetricTensor<rank,dim,Number>::n_independent_components ];
-      dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries );
-
-      SymmetricTensor<rank,dim,Number> global;
-      for (unsigned int i=0; i< n_entries; ++i)
-        global[ global.unrolled_to_component_indices(i) ] = global_entries[i];
-
-      return global;
-    }
-
-    template <typename T>
-    inline
-    T max (const T &t,
-           const MPI_Comm &mpi_communicator)
-    {
-      return internal::all_reduce(MPI_MAX, t, mpi_communicator);
-    }
-
-
     template <typename T, unsigned int N>
-    inline
     void max (const T (&values)[N],
               const MPI_Comm &mpi_communicator,
               T (&maxima)[N])
     {
-      internal::all_reduce(MPI_MAX, values, mpi_communicator, maxima);
-    }
-
-
-    template <typename T>
-    inline
-    void max (const std::vector<T> &values,
-              const MPI_Comm       &mpi_communicator,
-              std::vector<T>       &maxima)
-    {
-      internal::all_reduce(MPI_MAX, values, mpi_communicator, maxima);
-    }
-
-
-    template <typename T>
-    inline
-    T min (const T &t,
-           const MPI_Comm &mpi_communicator)
-    {
-      return internal::all_reduce(MPI_MIN, t, mpi_communicator);
+      internal::all_reduce(MPI_MAX, values, mpi_communicator, maxima, N);
     }
 
-
     template <typename T, unsigned int N>
-    inline
     void min (const T (&values)[N],
               const MPI_Comm &mpi_communicator,
               T (&minima)[N])
     {
-      internal::all_reduce(MPI_MIN, values, mpi_communicator, minima);
-    }
-
-
-    template <typename T>
-    inline
-    void min (const std::vector<T> &values,
-              const MPI_Comm       &mpi_communicator,
-              std::vector<T>       &minima)
-    {
-      internal::all_reduce(MPI_MIN, values, mpi_communicator, minima);
+      internal::all_reduce(MPI_MIN, values, mpi_communicator, minima, N);
     }
+#endif
 
 
     inline
diff --git a/include/deal.II/base/mpi.templates.h b/include/deal.II/base/mpi.templates.h
new file mode 100644 (file)
index 0000000..f4108f2
--- /dev/null
@@ -0,0 +1,263 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2011 - 2016 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.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii__mpi_templates_h
+#define dealii__mpi_templates_h
+
+#include <deal.II/base/config.h>
+#include <deal.II/base/exceptions.h>
+#include <deal.II/base/mpi.h>
+#include <deal.II/base/tensor.h>
+#include <deal.II/base/symmetric_tensor.h>
+#include <deal.II/lac/vector.h>
+
+#include <vector>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Utilities
+{
+  namespace MPI
+  {
+    namespace internal
+    {
+#ifdef DEAL_II_WITH_MPI
+      /**
+       * Return the corresponding MPI data type id for the argument given.
+       */
+      MPI_Datatype mpi_type_id (const int *)
+      {
+        return MPI_INT;
+      }
+
+
+      MPI_Datatype mpi_type_id (const long int *)
+      {
+        return MPI_LONG;
+      }
+
+
+      MPI_Datatype mpi_type_id (const unsigned int *)
+      {
+        return MPI_UNSIGNED;
+      }
+
+
+      MPI_Datatype mpi_type_id (const unsigned long int *)
+      {
+        return MPI_UNSIGNED_LONG;
+      }
+
+
+      MPI_Datatype mpi_type_id (const unsigned long long int *)
+      {
+        return MPI_UNSIGNED_LONG_LONG;
+      }
+
+
+      MPI_Datatype mpi_type_id (const float *)
+      {
+        return MPI_FLOAT;
+      }
+
+
+      MPI_Datatype mpi_type_id (const double *)
+      {
+        return MPI_DOUBLE;
+      }
+
+
+      MPI_Datatype mpi_type_id (const long double *)
+      {
+        return MPI_LONG_DOUBLE;
+      }
+#endif
+
+      template <typename T>
+      void all_reduce (const MPI_Op      &mpi_op,
+                       const T *const    values,
+                       const MPI_Comm    &mpi_communicator,
+                       T                 *output,
+                       const std::size_t  size)
+      {
+#ifdef DEAL_II_WITH_MPI
+        if (job_supports_mpi())
+          {
+            MPI_Allreduce (values != output
+                           ?
+                           // TODO This const_cast is only needed for older
+                           // (e.g., openMPI 1.6, released in 2012)
+                           // implementations of MPI-2. It is not needed as of
+                           // MPI-3 and we should remove it at some point in
+                           // the future.
+                           const_cast<void *>(static_cast<const void *>(values))
+                           :
+                           MPI_IN_PLACE,
+                           static_cast<void *>(output),
+                           static_cast<int>(size),
+                           internal::mpi_type_id(values),
+                           mpi_op,
+                           mpi_communicator);
+          }
+        else
+#endif
+          {
+            (void)mpi_op;
+            (void)mpi_communicator;
+            for (std::size_t i=0; i<size; ++i)
+              output[i] = values[i];
+          }
+      }
+
+      template <typename T>
+      T all_reduce (const MPI_Op &mpi_op,
+                    const T &t,
+                    const MPI_Comm &mpi_communicator)
+      {
+        T output;
+        all_reduce(mpi_op, &t, mpi_communicator, &output, 1);
+        return output;
+      }
+
+      template <typename T>
+      void all_reduce (const MPI_Op &mpi_op,
+                       const std::vector<T> &values,
+                       const MPI_Comm       &mpi_communicator,
+                       std::vector<T>       &output)
+      {
+        Assert(values.size() == output.size(),
+               ExcDimensionMismatch(values.size(), output.size()));
+        all_reduce(mpi_op, &values[0], mpi_communicator, &output[0], values.size());
+      }
+
+      template <typename T>
+      void all_reduce (const MPI_Op    &mpi_op,
+                       const Vector<T> &values,
+                       const MPI_Comm  &mpi_communicator,
+                       Vector<T>  &output)
+      {
+        Assert(values.size() == output.size(),
+               ExcDimensionMismatch(values.size(), output.size()));
+        all_reduce(mpi_op, values.begin(), mpi_communicator, output.begin(), values.size());
+      }
+    }
+
+
+    template <typename T>
+    T sum (const T &t,
+           const MPI_Comm &mpi_communicator)
+    {
+      return internal::all_reduce(MPI_SUM, t, mpi_communicator);
+    }
+
+
+    template <typename T>
+    void sum (const std::vector<T> &values,
+              const MPI_Comm       &mpi_communicator,
+              std::vector<T>       &sums)
+    {
+      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums);
+    }
+
+    template <typename T>
+    void sum (const Vector<T> &values,
+              const MPI_Comm &mpi_communicator,
+              Vector<T> &sums)
+    {
+      internal::all_reduce(MPI_SUM, values, mpi_communicator, sums);
+    }
+
+
+    template <int rank, int dim, typename Number>
+    Tensor<rank,dim,Number>
+    sum (const Tensor<rank,dim,Number> &local,
+         const MPI_Comm &mpi_communicator)
+    {
+      const unsigned int n_entries = Tensor<rank,dim,Number>::n_independent_components;
+      Number entries[ Tensor<rank,dim,Number>::n_independent_components ];
+
+      for (unsigned int i=0; i< n_entries; ++i)
+        entries[i] = local[ local.unrolled_to_component_indices(i) ];
+
+      Number global_entries[ Tensor<rank,dim,Number>::n_independent_components ];
+      dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries );
+
+      Tensor<rank,dim,Number> global;
+      for (unsigned int i=0; i< n_entries; ++i)
+        global[ global.unrolled_to_component_indices(i) ] = global_entries[i];
+
+      return global;
+    }
+
+    template <int rank, int dim, typename Number>
+    SymmetricTensor<rank,dim,Number>
+    sum (const SymmetricTensor<rank,dim,Number> &local,
+         const MPI_Comm &mpi_communicator)
+    {
+      const unsigned int n_entries = SymmetricTensor<rank,dim,Number>::n_independent_components;
+      Number entries[ SymmetricTensor<rank,dim,Number>::n_independent_components ];
+
+      for (unsigned int i=0; i< n_entries; ++i)
+        entries[i] = local[ local.unrolled_to_component_indices(i) ];
+
+      Number global_entries[ SymmetricTensor<rank,dim,Number>::n_independent_components ];
+      dealii::Utilities::MPI::sum( entries, mpi_communicator, global_entries );
+
+      SymmetricTensor<rank,dim,Number> global;
+      for (unsigned int i=0; i< n_entries; ++i)
+        global[ global.unrolled_to_component_indices(i) ] = global_entries[i];
+
+      return global;
+    }
+
+    template <typename T>
+    T max (const T &t,
+           const MPI_Comm &mpi_communicator)
+    {
+      return internal::all_reduce(MPI_MAX, t, mpi_communicator);
+    }
+
+
+    template <typename T>
+    void max (const std::vector<T> &values,
+              const MPI_Comm       &mpi_communicator,
+              std::vector<T>       &maxima)
+    {
+      internal::all_reduce(MPI_MAX, values, mpi_communicator, maxima);
+    }
+
+
+    template <typename T>
+    T min (const T &t,
+           const MPI_Comm &mpi_communicator)
+    {
+      return internal::all_reduce(MPI_MIN, t, mpi_communicator);
+    }
+
+
+    template <typename T>
+    void min (const std::vector<T> &values,
+              const MPI_Comm       &mpi_communicator,
+              std::vector<T>       &minima)
+    {
+      internal::all_reduce(MPI_MIN, values, mpi_communicator, minima);
+    }
+  } // end of namespace MPI
+} // end of namespace Utilities
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 587b1dada98c496c88ae2f409c281076c9ed4be1..9df7b72fb6ba7da9e1f6e09dd3f15a74e758abc2 100644 (file)
@@ -75,6 +75,7 @@ SET(_inst
   data_out_base.inst.in
   function.inst.in
   function_time.inst.in
+  mpi.inst.in
   polynomials_rannacher_turek.inst.in
   tensor_function.inst.in
   time_stepping.inst.in
index 858a621f4f3d65b2be351b2b20a8ca257f5c12ed..4c86df4589461479a7b21d8e1cbda835ce238783 100644 (file)
@@ -15,6 +15,7 @@
 
 
 #include <deal.II/base/mpi.h>
+#include <deal.II/base/mpi.templates.h>
 #include <deal.II/base/utilities.h>
 #include <deal.II/base/exceptions.h>
 #include <deal.II/lac/vector_memory.h>
@@ -506,10 +507,9 @@ namespace Utilities
         }
 #endif
     }
-
-
+#include "mpi.inst"
+    }
   } // end of namespace MPI
-
 } // end of namespace Utilities
 
 DEAL_II_NAMESPACE_CLOSE
diff --git a/source/base/mpi.inst.in b/source/base/mpi.inst.in
new file mode 100644 (file)
index 0000000..531f86b
--- /dev/null
@@ -0,0 +1,60 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 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.
+//
+// ---------------------------------------------------------------------
+
+
+for (S : MPI_SCALARS)
+{
+  template
+  void sum<S> (const Vector<S> &, const MPI_Comm &, Vector<S> &);
+
+  template
+  S sum<S> (const S &, const MPI_Comm &);
+
+  template
+  void sum<S> (const std::vector<S> &, const MPI_Comm &, std::vector<S> &);
+
+  template
+  S max<S> (const S &, const MPI_Comm &);
+
+  template
+  void max<S> (const std::vector<S> &, const MPI_Comm &, std::vector<S> &);
+
+  template
+  S min<S> (const S &, const MPI_Comm &);
+
+  template
+  void min<S> (const std::vector<S> &, const MPI_Comm &, std::vector<S> &);
+}
+
+
+for (S : REAL_SCALARS; rank: RANKS; dim : SPACE_DIMENSIONS)
+{
+  template
+  Tensor<rank,dim,S>
+  sum<rank,dim,S> (const Tensor<rank,dim,S> &, const MPI_Comm &);
+}
+
+
+// Recall that SymmetricTensor only makes sense for rank 2 and 4
+for (S : REAL_SCALARS; dim : SPACE_DIMENSIONS)
+{
+  template
+  SymmetricTensor<2,dim,S>
+  sum<2,dim,S> (const SymmetricTensor<2,dim,S> &, const MPI_Comm &);
+
+  template
+  SymmetricTensor<4,dim,S>
+  sum<4,dim,S> (const SymmetricTensor<4,dim,S> &, const MPI_Comm &);
+}
index c72a61c3eea64dc17fb0da29c1b57dd4adbe5aff..9fab96ae2b45aea8632adf69168a699da529fa2d 100644 (file)
@@ -22,6 +22,7 @@
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
 #include <deal.II/base/mpi.h>
+#include <deal.II/base/mpi.templates.h>
 #include <deal.II/distributed/tria.h>
 #include <deal.II/grid/grid_generator.h>
 #include <deal.II/dofs/dof_handler.h>
index d0ecc31957f268e7394eba873fe73178fef50cea..99de134f94dff3e4901fd9081a74489a82427e19 100644 (file)
@@ -22,6 +22,7 @@
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
 #include <deal.II/base/mpi.h>
+#include <deal.II/base/mpi.templates.h>
 #include <deal.II/distributed/tria.h>
 #include <deal.II/grid/grid_generator.h>
 #include <deal.II/dofs/dof_handler.h>
index 7a1dd16a803f385230220dc2efc40fd65e2c2945..d8ae92099055592c148d43b2515f4708beb6b4a4 100644 (file)
@@ -21,6 +21,7 @@
 #include <deal.II/base/config.h>
 #include <deal.II/base/job_identifier.h>
 #include <deal.II/base/logstream.h>
+#include <deal.II/base/mpi.h>
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/utilities.h>
 #include <deal.II/base/thread_management.h>

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.