Reference documentation for deal.II version 8.5.1
mpi.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__mpi_h
17 #define dealii__mpi_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <vector>
22 
23 #if !defined(DEAL_II_WITH_MPI) && !defined(DEAL_II_WITH_PETSC)
24 // without MPI, we would still like to use
25 // some constructs with MPI data
26 // types. Therefore, create some dummies
27 typedef int MPI_Comm;
28 const int MPI_COMM_SELF = 0;
29 typedef int MPI_Datatype;
30 typedef int MPI_Op;
31 
32 static const int MPI_MIN = 0;
33 static const int MPI_MAX = 0;
34 static const int MPI_SUM = 0;
35 #endif
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 //Forward type declarations to allow MPI sums over tensorial types
41 template <int rank, int dim, typename Number> class Tensor;
42 template <int rank, int dim, typename Number> class SymmetricTensor;
43 
44 //Forward type declaration to allow MPI sums over Vector<number> type
45 template <typename Number> class Vector;
46 
47 
48 namespace Utilities
49 {
57  namespace MPI
58  {
64  unsigned int n_mpi_processes (const MPI_Comm &mpi_communicator);
65 
74  unsigned int this_mpi_process (const MPI_Comm &mpi_communicator);
75 
97  std::vector<unsigned int>
98  compute_point_to_point_communication_pattern (const MPI_Comm &mpi_comm,
99  const std::vector<unsigned int> &destinations);
100 
114  MPI_Comm duplicate_communicator (const MPI_Comm &mpi_communicator);
115 
135  template <typename T>
136  T sum (const T &t,
137  const MPI_Comm &mpi_communicator);
138 
147  template <typename T, unsigned int N>
148  void sum (const T (&values)[N],
149  const MPI_Comm &mpi_communicator,
150  T (&sums)[N]);
151 
159  template <typename T>
160  void sum (const std::vector<T> &values,
161  const MPI_Comm &mpi_communicator,
162  std::vector<T> &sums);
163 
170  template <typename T>
171  void sum (const Vector<T> &values,
172  const MPI_Comm &mpi_communicator,
173  Vector<T> &sums);
174 
175 
181  template <int rank, int dim, typename Number>
184  const MPI_Comm &mpi_communicator);
185 
191  template <int rank, int dim, typename Number>
193  sum (const Tensor<rank,dim,Number> &local,
194  const MPI_Comm &mpi_communicator);
195 
215  template <typename T>
216  T max (const T &t,
217  const MPI_Comm &mpi_communicator);
218 
227  template <typename T, unsigned int N>
228  void max (const T (&values)[N],
229  const MPI_Comm &mpi_communicator,
230  T (&maxima)[N]);
231 
240  template <typename T>
241  void max (const std::vector<T> &values,
242  const MPI_Comm &mpi_communicator,
243  std::vector<T> &maxima);
244 
264  template <typename T>
265  T min (const T &t,
266  const MPI_Comm &mpi_communicator);
267 
276  template <typename T, unsigned int N>
277  void min (const T (&values)[N],
278  const MPI_Comm &mpi_communicator,
279  T (&minima)[N]);
280 
289  template <typename T>
290  void min (const std::vector<T> &values,
291  const MPI_Comm &mpi_communicator,
292  std::vector<T> &minima);
293 
294 
298  struct MinMaxAvg
299  {
300  // Note: We assume a POD property of this struct in the MPI calls in
301  // min_max_avg
302  double sum;
303  double min;
304  double max;
305  unsigned int min_index;
306  unsigned int max_index;
307  double avg;
308  };
309 
324  MinMaxAvg
325  min_max_avg (const double my_value,
326  const MPI_Comm &mpi_communicator);
327 
328 
329 
353  {
354  public:
400  MPI_InitFinalize (int &argc,
401  char ** &argv,
402  const unsigned int max_num_threads = numbers::invalid_unsigned_int);
403 
409  };
410 
422  bool job_supports_mpi ();
423 
424 #ifndef DOXYGEN
425  // declaration for an internal function that lives in mpi.templates.h
426  namespace internal
427  {
428  template <typename T>
429  void all_reduce (const MPI_Op &mpi_op,
430  const T *const values,
431  const MPI_Comm &mpi_communicator,
432  T *output,
433  const std::size_t size);
434  }
435 
436  // Since these depend on N they must live in the header file
437  template <typename T, unsigned int N>
438  void sum (const T (&values)[N],
439  const MPI_Comm &mpi_communicator,
440  T (&sums)[N])
441  {
442  internal::all_reduce(MPI_SUM, values, mpi_communicator, sums, N);
443  }
444 
445  template <typename T, unsigned int N>
446  void max (const T (&values)[N],
447  const MPI_Comm &mpi_communicator,
448  T (&maxima)[N])
449  {
450  internal::all_reduce(MPI_MAX, values, mpi_communicator, maxima, N);
451  }
452 
453  template <typename T, unsigned int N>
454  void min (const T (&values)[N],
455  const MPI_Comm &mpi_communicator,
456  T (&minima)[N])
457  {
458  internal::all_reduce(MPI_MIN, values, mpi_communicator, minima, N);
459  }
460 #endif
461  } // end of namespace MPI
462 } // end of namespace Utilities
463 
464 
465 DEAL_II_NAMESPACE_CLOSE
466 
467 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:170
MPI_InitFinalize(int &argc, char **&argv, const unsigned int max_num_threads=numbers::invalid_unsigned_int)
Definition: mpi.cc:297
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:63
Definition: mpi.h:48
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:83
Definition: mpi.h:41
T min(const T &t, const MPI_Comm &mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:73
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:194
std::vector< unsigned int > compute_point_to_point_communication_pattern(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
Definition: mpi.cc:93
bool job_supports_mpi()
Definition: mpi.cc:519
T max(const T &t, const MPI_Comm &mpi_communicator)