]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Also allow MPI::sum on complex numbers
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sat, 11 Jun 2016 20:23:41 +0000 (22:23 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Tue, 14 Jun 2016 15:10:11 +0000 (17:10 +0200)
include/deal.II/base/mpi.templates.h
source/base/mpi.inst.in

index f4108f21420dc4ae0016d89af418f893ac3e4eff..a5806cc217b1c9639be0b610daf14db60bb648fa 100644 (file)
@@ -121,6 +121,43 @@ namespace Utilities
           }
       }
 
+      template <typename T>
+      void all_reduce (const MPI_Op                   &mpi_op,
+                       const std::complex<T> *const    values,
+                       const MPI_Comm                 &mpi_communicator,
+                       std::complex<T>                *output,
+                       const std::size_t               size)
+      {
+#ifdef DEAL_II_WITH_MPI
+        if (job_supports_mpi())
+          {
+            T dummy_selector;
+            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*2),
+                           internal::mpi_type_id(&dummy_selector),
+                           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,
index 531f86beb14c0fe2764355d60cf481db00ef2d40..19a06a032828ab7826357554f29a2c27064a7daf 100644 (file)
@@ -39,6 +39,19 @@ for (S : MPI_SCALARS)
 }
 
 
+for (S : COMPLEX_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> &);
+}
+
+
 for (S : REAL_SCALARS; rank: RANKS; dim : SPACE_DIMENSIONS)
 {
   template

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.