Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mpi.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2020 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mpi_h
17 #define dealii_mpi_h
18 
19 #include <deal.II/base/config.h>
20 
22 #include <deal.II/base/mpi_tags.h>
23 #include <deal.II/base/numbers.h>
24 
25 #include <map>
26 #include <numeric>
27 #include <set>
28 #include <vector>
29 
30 #if !defined(DEAL_II_WITH_MPI) && !defined(DEAL_II_WITH_PETSC)
31 // without MPI, we would still like to use
32 // some constructs with MPI data
33 // types. Therefore, create some dummies
34 using MPI_Comm = int;
35 using MPI_Request = int;
36 using MPI_Datatype = int;
37 using MPI_Op = int;
38 # ifndef MPI_COMM_WORLD
39 # define MPI_COMM_WORLD 0
40 # endif
41 # ifndef MPI_COMM_SELF
42 # define MPI_COMM_SELF 0
43 # endif
44 # ifndef MPI_REQUEST_NULL
45 # define MPI_REQUEST_NULL 0
46 # endif
47 # ifndef MPI_MIN
48 # define MPI_MIN 0
49 # endif
50 # ifndef MPI_MAX
51 # define MPI_MAX 0
52 # endif
53 # ifndef MPI_SUM
54 # define MPI_SUM 0
55 # endif
56 #endif
57 
58 
59 
73 #ifdef DEAL_II_WITH_MPI
74 # if DEAL_II_MPI_VERSION_GTE(3, 0)
75 
76 # define DEAL_II_MPI_CONST_CAST(expr) (expr)
77 
78 # else
79 
80 # include <type_traits>
81 
82 # define DEAL_II_MPI_CONST_CAST(expr) \
83  const_cast<typename std::remove_const< \
84  typename std::remove_pointer<decltype(expr)>::type>::type *>(expr)
85 
86 # endif
87 #endif
88 
89 
90 
92 
93 
94 // Forward type declarations to allow MPI sums over tensorial types
95 #ifndef DOXYGEN
96 template <int rank, int dim, typename Number>
97 class Tensor;
98 template <int rank, int dim, typename Number>
99 class SymmetricTensor;
100 template <typename Number>
101 class SparseMatrix;
102 class IndexSet;
103 #endif
104 
105 namespace Utilities
106 {
119  IndexSet
120  create_evenly_distributed_partitioning(const unsigned int my_partition_id,
121  const unsigned int n_partitions,
122  const IndexSet::size_type total_size);
123 
131  namespace MPI
132  {
141  unsigned int
142  n_mpi_processes(const MPI_Comm &mpi_communicator);
143 
152  unsigned int
153  this_mpi_process(const MPI_Comm &mpi_communicator);
154 
176  std::vector<unsigned int>
178  const MPI_Comm & mpi_comm,
179  const std::vector<unsigned int> &destinations);
180 
200  unsigned int
202  const MPI_Comm & mpi_comm,
203  const std::vector<unsigned int> &destinations);
204 
221  MPI_Comm
222  duplicate_communicator(const MPI_Comm &mpi_communicator);
223 
233  void
234  free_communicator(MPI_Comm &mpi_communicator);
235 
249  {
250  public:
254  explicit DuplicatedCommunicator(const MPI_Comm &communicator)
255  : comm(duplicate_communicator(communicator))
256  {}
257 
262 
267  {
269  }
270 
274  const MPI_Comm &operator*() const
275  {
276  return comm;
277  }
278 
279 
284  operator=(const DuplicatedCommunicator &) = delete;
285 
286  private:
291  };
292 
323  {
324  public:
331  {
332  public:
337  : mutex(mutex)
338  , comm(comm)
339  {
340  mutex.lock(comm);
341  }
342 
347  {
348  mutex.unlock(comm);
349  }
350 
351  private:
359  const MPI_Comm comm;
360  };
361 
365  explicit CollectiveMutex();
366 
371 
378  void
379  lock(MPI_Comm comm);
380 
387  void
388  unlock(MPI_Comm comm);
389 
390  private:
394  bool locked;
395 
400  };
401 
402 
403 
431 #ifdef DEAL_II_WITH_MPI
432  int
433  create_group(const MPI_Comm & comm,
434  const MPI_Group &group,
435  const int tag,
436  MPI_Comm * new_comm);
437 #endif
438 
447  std::vector<IndexSet>
449  const IndexSet::size_type local_size);
450 
458  IndexSet
460  const MPI_Comm & comm,
461  const IndexSet::size_type total_size);
462 
463 #ifdef DEAL_II_WITH_MPI
464 
479  template <class Iterator, typename Number = long double>
480  std::pair<Number, typename numbers::NumberTraits<Number>::real_type>
482  const Iterator end,
483  const MPI_Comm &comm);
484 #endif
485 
505  template <typename T>
506  T
507  sum(const T &t, const MPI_Comm &mpi_communicator);
508 
518  template <typename T, typename U>
519  void
520  sum(const T &values, const MPI_Comm &mpi_communicator, U &sums);
521 
531  template <typename T>
532  void
533  sum(const ArrayView<const T> &values,
534  const MPI_Comm & mpi_communicator,
535  const ArrayView<T> & sums);
536 
542  template <int rank, int dim, typename Number>
545  const MPI_Comm & mpi_communicator);
546 
552  template <int rank, int dim, typename Number>
554  sum(const Tensor<rank, dim, Number> &local,
555  const MPI_Comm & mpi_communicator);
556 
565  template <typename Number>
566  void
567  sum(const SparseMatrix<Number> &local,
568  const MPI_Comm & mpi_communicator,
569  SparseMatrix<Number> & global);
570 
590  template <typename T>
591  T
592  max(const T &t, const MPI_Comm &mpi_communicator);
593 
603  template <typename T, typename U>
604  void
605  max(const T &values, const MPI_Comm &mpi_communicator, U &maxima);
606 
616  template <typename T>
617  void
618  max(const ArrayView<const T> &values,
619  const MPI_Comm & mpi_communicator,
620  const ArrayView<T> & maxima);
621 
641  template <typename T>
642  T
643  min(const T &t, const MPI_Comm &mpi_communicator);
644 
654  template <typename T, typename U>
655  void
656  min(const T &values, const MPI_Comm &mpi_communicator, U &minima);
657 
667  template <typename T>
668  void
669  min(const ArrayView<const T> &values,
670  const MPI_Comm & mpi_communicator,
671  const ArrayView<T> & minima);
672 
687  struct MinMaxAvg
688  {
693  double sum;
694 
699  double min;
700 
705  double max;
706 
715  unsigned int min_index;
716 
725  unsigned int max_index;
726 
731  double avg;
732  };
733 
748  MinMaxAvg
749  min_max_avg(const double my_value, const MPI_Comm &mpi_communicator);
750 
762  std::vector<MinMaxAvg>
763  min_max_avg(const std::vector<double> &my_value,
764  const MPI_Comm & mpi_communicator);
765 
766 
779  void
780  min_max_avg(const ArrayView<const double> &my_values,
781  const ArrayView<MinMaxAvg> & result,
782  const MPI_Comm & mpi_communicator);
783 
784 
829  {
830  public:
877  int & argc,
878  char **& argv,
879  const unsigned int max_num_threads = numbers::invalid_unsigned_int);
880 
886 
913  static void
914  register_request(MPI_Request &request);
915 
919  static void
921 
922  private:
926  static std::set<MPI_Request *> requests;
927  };
928 
940  bool
942 
962  template <typename T>
963  std::map<unsigned int, T>
964  some_to_some(const MPI_Comm & comm,
965  const std::map<unsigned int, T> &objects_to_send);
966 
982  template <typename T>
983  std::vector<T>
984  all_gather(const MPI_Comm &comm, const T &object_to_send);
985 
1003  template <typename T>
1004  std::vector<T>
1005  gather(const MPI_Comm & comm,
1006  const T & object_to_send,
1007  const unsigned int root_process = 0);
1008 
1053  std::vector<unsigned int>
1054  compute_index_owner(const IndexSet &owned_indices,
1055  const IndexSet &indices_to_look_up,
1056  const MPI_Comm &comm);
1057 
1065  template <typename T>
1066  std::vector<T>
1067  compute_set_union(const std::vector<T> &vec, const MPI_Comm &comm);
1068 
1072  template <typename T>
1073  std::set<T>
1074  compute_set_union(const std::set<T> &set, const MPI_Comm &comm);
1075 
1076 #ifndef DOXYGEN
1077  // declaration for an internal function that lives in mpi.templates.h
1078  namespace internal
1079  {
1080  template <typename T>
1081  void
1082  all_reduce(const MPI_Op & mpi_op,
1083  const ArrayView<const T> &values,
1084  const MPI_Comm & mpi_communicator,
1085  const ArrayView<T> & output);
1086  }
1087 
1088  // Since these depend on N they must live in the header file
1089  template <typename T, unsigned int N>
1090  void
1091  sum(const T (&values)[N], const MPI_Comm &mpi_communicator, T (&sums)[N])
1092  {
1093  internal::all_reduce(MPI_SUM,
1094  ArrayView<const T>(values, N),
1095  mpi_communicator,
1096  ArrayView<T>(sums, N));
1097  }
1098 
1099  template <typename T, unsigned int N>
1100  void
1101  max(const T (&values)[N], const MPI_Comm &mpi_communicator, T (&maxima)[N])
1102  {
1103  internal::all_reduce(MPI_MAX,
1104  ArrayView<const T>(values, N),
1105  mpi_communicator,
1106  ArrayView<T>(maxima, N));
1107  }
1108 
1109  template <typename T, unsigned int N>
1110  void
1111  min(const T (&values)[N], const MPI_Comm &mpi_communicator, T (&minima)[N])
1112  {
1113  internal::all_reduce(MPI_MIN,
1114  ArrayView<const T>(values, N),
1115  mpi_communicator,
1116  ArrayView<T>(minima, N));
1117  }
1118 
1119  template <typename T>
1120  std::map<unsigned int, T>
1121  some_to_some(const MPI_Comm & comm,
1122  const std::map<unsigned int, T> &objects_to_send)
1123  {
1124 # ifndef DEAL_II_WITH_MPI
1125  (void)comm;
1126  Assert(objects_to_send.size() < 2,
1127  ExcMessage("Cannot send to more than one processor."));
1128  Assert(objects_to_send.find(0) != objects_to_send.end() ||
1129  objects_to_send.size() == 0,
1130  ExcMessage("Can only send to myself or to nobody."));
1131  return objects_to_send;
1132 # else
1133  const auto my_proc = this_mpi_process(comm);
1134 
1135  std::map<unsigned int, T> received_objects;
1136 
1137  std::vector<unsigned int> send_to;
1138  send_to.reserve(objects_to_send.size());
1139  for (const auto &m : objects_to_send)
1140  if (m.first == my_proc)
1141  received_objects[my_proc] = m.second;
1142  else
1143  send_to.emplace_back(m.first);
1144 
1145  const unsigned int n_point_point_communications =
1147 
1148  // Protect the following communication:
1149  static CollectiveMutex mutex;
1150  CollectiveMutex::ScopedLock lock(mutex, comm);
1151 
1152  // If we have something to send, or we expect something from other
1153  // processors, we need to visit one of the two scopes below. Otherwise,
1154  // no other action is required by this mpi process, and we can safely
1155  // return.
1156  if (send_to.size() == 0 && n_point_point_communications == 0)
1157  return received_objects;
1158 
1159  const int mpi_tag =
1161 
1162  // Sending buffers
1163  std::vector<std::vector<char>> buffers_to_send(send_to.size());
1164  std::vector<MPI_Request> buffer_send_requests(send_to.size());
1165  {
1166  unsigned int i = 0;
1167  for (const auto &rank_obj : objects_to_send)
1168  if (rank_obj.first != my_proc)
1169  {
1170  const auto &rank = rank_obj.first;
1171  buffers_to_send[i] = Utilities::pack(rank_obj.second);
1172  const int ierr = MPI_Isend(buffers_to_send[i].data(),
1173  buffers_to_send[i].size(),
1174  MPI_CHAR,
1175  rank,
1176  mpi_tag,
1177  comm,
1178  &buffer_send_requests[i]);
1179  AssertThrowMPI(ierr);
1180  ++i;
1181  }
1182  }
1183 
1184  // Fill the output map
1185  {
1186  std::vector<char> buffer;
1187  // We do this on a first come/first served basis
1188  for (unsigned int i = 0; i < n_point_point_communications; ++i)
1189  {
1190  // Probe what's going on. Take data from the first available sender
1191  MPI_Status status;
1192  int ierr = MPI_Probe(MPI_ANY_SOURCE, mpi_tag, comm, &status);
1193  AssertThrowMPI(ierr);
1194 
1195  // Length of the message
1196  int len;
1197  ierr = MPI_Get_count(&status, MPI_CHAR, &len);
1198  AssertThrowMPI(ierr);
1199  buffer.resize(len);
1200 
1201  // Source rank
1202  const unsigned int rank = status.MPI_SOURCE;
1203 
1204  // Actually receive the message
1205  ierr = MPI_Recv(buffer.data(),
1206  len,
1207  MPI_CHAR,
1208  status.MPI_SOURCE,
1209  status.MPI_TAG,
1210  comm,
1211  MPI_STATUS_IGNORE);
1212  AssertThrowMPI(ierr);
1213  Assert(received_objects.find(rank) == received_objects.end(),
1215  "I should not receive again from this rank"));
1216  received_objects[rank] = Utilities::unpack<T>(buffer);
1217  }
1218  }
1219 
1220  // Wait to have sent all objects.
1221  const int ierr = MPI_Waitall(send_to.size(),
1222  buffer_send_requests.data(),
1223  MPI_STATUSES_IGNORE);
1224  AssertThrowMPI(ierr);
1225 
1226  return received_objects;
1227 # endif // deal.II with MPI
1228  }
1229 
1230  template <typename T>
1231  std::vector<T>
1232  all_gather(const MPI_Comm &comm, const T &object)
1233  {
1234  if (job_supports_mpi() == false)
1235  return {object};
1236 
1237 # ifndef DEAL_II_WITH_MPI
1238  (void)comm;
1239  std::vector<T> v(1, object);
1240  return v;
1241 # else
1242  const auto n_procs = ::Utilities::MPI::n_mpi_processes(comm);
1243 
1244  std::vector<char> buffer = Utilities::pack(object);
1245 
1246  int n_local_data = buffer.size();
1247 
1248  // Vector to store the size of loc_data_array for every process
1249  std::vector<int> size_all_data(n_procs, 0);
1250 
1251  // Exchanging the size of each buffer
1252  MPI_Allgather(
1253  &n_local_data, 1, MPI_INT, size_all_data.data(), 1, MPI_INT, comm);
1254 
1255  // Now computing the displacement, relative to recvbuf,
1256  // at which to store the incoming buffer
1257  std::vector<int> rdispls(n_procs);
1258  rdispls[0] = 0;
1259  for (unsigned int i = 1; i < n_procs; ++i)
1260  rdispls[i] = rdispls[i - 1] + size_all_data[i - 1];
1261 
1262  // Step 3: exchange the buffer:
1263  std::vector<char> received_unrolled_buffer(rdispls.back() +
1264  size_all_data.back());
1265 
1266  MPI_Allgatherv(buffer.data(),
1267  n_local_data,
1268  MPI_CHAR,
1269  received_unrolled_buffer.data(),
1270  size_all_data.data(),
1271  rdispls.data(),
1272  MPI_CHAR,
1273  comm);
1274 
1275  std::vector<T> received_objects(n_procs);
1276  for (unsigned int i = 0; i < n_procs; ++i)
1277  {
1278  std::vector<char> local_buffer(received_unrolled_buffer.begin() +
1279  rdispls[i],
1280  received_unrolled_buffer.begin() +
1281  rdispls[i] + size_all_data[i]);
1282  received_objects[i] = Utilities::unpack<T>(local_buffer);
1283  }
1284 
1285  return received_objects;
1286 # endif
1287  }
1288 
1289  template <typename T>
1290  std::vector<T>
1291  gather(const MPI_Comm & comm,
1292  const T & object_to_send,
1293  const unsigned int root_process)
1294  {
1295 # ifndef DEAL_II_WITH_MPI
1296  (void)comm;
1297  (void)root_process;
1298  std::vector<T> v(1, object_to_send);
1299  return v;
1300 # else
1301  const auto n_procs = ::Utilities::MPI::n_mpi_processes(comm);
1302  const auto my_rank = ::Utilities::MPI::this_mpi_process(comm);
1303 
1304  AssertIndexRange(root_process, n_procs);
1305 
1306  std::vector<char> buffer = Utilities::pack(object_to_send);
1307  int n_local_data = buffer.size();
1308 
1309  // Vector to store the size of loc_data_array for every process
1310  // only the root process needs to allocate memory for that purpose
1311  std::vector<int> size_all_data;
1312  if (my_rank == root_process)
1313  size_all_data.resize(n_procs, 0);
1314 
1315  // Exchanging the size of each buffer
1316  int ierr = MPI_Gather(&n_local_data,
1317  1,
1318  MPI_INT,
1319  size_all_data.data(),
1320  1,
1321  MPI_INT,
1322  root_process,
1323  comm);
1324  AssertThrowMPI(ierr);
1325 
1326  // Now computing the displacement, relative to recvbuf,
1327  // at which to store the incoming buffer; only for root
1328  std::vector<int> rdispls;
1329  if (my_rank == root_process)
1330  {
1331  rdispls.resize(n_procs, 0);
1332  for (unsigned int i = 1; i < n_procs; ++i)
1333  rdispls[i] = rdispls[i - 1] + size_all_data[i - 1];
1334  }
1335  // exchange the buffer:
1336  std::vector<char> received_unrolled_buffer;
1337  if (my_rank == root_process)
1338  received_unrolled_buffer.resize(rdispls.back() + size_all_data.back());
1339 
1340  ierr = MPI_Gatherv(buffer.data(),
1341  n_local_data,
1342  MPI_CHAR,
1343  received_unrolled_buffer.data(),
1344  size_all_data.data(),
1345  rdispls.data(),
1346  MPI_CHAR,
1347  root_process,
1348  comm);
1349  AssertThrowMPI(ierr);
1350 
1351  std::vector<T> received_objects;
1352 
1353  if (my_rank == root_process)
1354  {
1355  received_objects.resize(n_procs);
1356 
1357  for (unsigned int i = 0; i < n_procs; ++i)
1358  {
1359  const std::vector<char> local_buffer(
1360  received_unrolled_buffer.begin() + rdispls[i],
1361  received_unrolled_buffer.begin() + rdispls[i] +
1362  size_all_data[i]);
1363  received_objects[i] = Utilities::unpack<T>(local_buffer);
1364  }
1365  }
1366  return received_objects;
1367 # endif
1368  }
1369 
1370 
1371 # ifdef DEAL_II_WITH_MPI
1372  template <class Iterator, typename Number>
1373  std::pair<Number, typename numbers::NumberTraits<Number>::real_type>
1375  const Iterator end,
1376  const MPI_Comm &comm)
1377  {
1378  // below we do simple and straight-forward implementation. More elaborate
1379  // options are:
1380  // http://dx.doi.org/10.1145/2807591.2807644 section 3.1.2
1381  // https://en.wikipedia.org/wiki/Algorithms_for_calculating_variance#Welford's_online_algorithm
1382  // https://en.wikipedia.org/wiki/Algorithms_for_calculating_variance#Online
1383  using Std = typename numbers::NumberTraits<Number>::real_type;
1384  const Number sum = std::accumulate(begin, end, Number(0.));
1385 
1386  const auto size = Utilities::MPI::sum(std::distance(begin, end), comm);
1387  Assert(size > 0, ExcDivideByZero());
1388  const Number mean =
1389  Utilities::MPI::sum(sum, comm) / static_cast<Std>(size);
1390  Std sq_sum = 0.;
1391  std::for_each(begin, end, [&mean, &sq_sum](const Number &v) {
1393  });
1394  sq_sum = Utilities::MPI::sum(sq_sum, comm);
1395  return std::make_pair(mean,
1396  std::sqrt(sq_sum / static_cast<Std>(size - 1)));
1397  }
1398 # endif
1399 
1400 #endif
1401  } // end of namespace MPI
1402 } // end of namespace Utilities
1403 
1404 
1406 
1407 #endif
array_view.h
IndexSet
Definition: index_set.h:74
Utilities::MPI::MinMaxAvg::sum
double sum
Definition: mpi.h:693
SymmetricTensor
Definition: symmetric_tensor.h:611
Utilities::MPI::create_evenly_distributed_partitioning
IndexSet create_evenly_distributed_partitioning(const MPI_Comm &comm, const IndexSet::size_type total_size)
Definition: mpi.cc:264
Utilities::MPI::DuplicatedCommunicator::operator=
DuplicatedCommunicator & operator=(const DuplicatedCommunicator &)=delete
Utilities::MPI::MinMaxAvg::min_index
unsigned int min_index
Definition: mpi.h:715
ArrayView
Definition: array_view.h:77
SparseMatrix
Definition: sparse_matrix.h:497
LAPACKSupport::U
static const char U
Definition: lapack_support.h:167
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
numbers::NumberTraits::abs_square
static constexpr std::enable_if< std::is_same< Dummy, number >::value &&is_cuda_compatible< Dummy >::value, real_type >::type abs_square(const number &x)
Definition: numbers.h:587
Utilities::MPI::DuplicatedCommunicator
Definition: mpi.h:248
Utilities::create_evenly_distributed_partitioning
IndexSet create_evenly_distributed_partitioning(const unsigned int my_partition_id, const unsigned int n_partitions, const IndexSet::size_type total_size)
Definition: mpi.cc:71
Utilities::MPI::MPI_InitFinalize::~MPI_InitFinalize
~MPI_InitFinalize()
Definition: mpi.cc:927
MPI_Comm
StandardExceptions::ExcDivideByZero
static ::ExceptionBase & ExcDivideByZero()
Utilities::pack
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
Definition: utilities.h:1209
Utilities::MPI::CollectiveMutex::ScopedLock::mutex
CollectiveMutex & mutex
Definition: mpi.h:355
Utilities::MPI::compute_n_point_to_point_communications
unsigned int compute_n_point_to_point_communications(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
Definition: mpi.cc:479
Utilities::MPI::internal::Tags::compute_point_to_point_communication_pattern
@ compute_point_to_point_communication_pattern
Utilities::MPI::compute_point_to_point_communication_pattern()
Definition: mpi_tags.h:57
LinearAlgebra::CUDAWrappers::kernel::set
__global__ void set(Number *val, const Number s, const size_type N)
Utilities::MPI::MinMaxAvg::max
double max
Definition: mpi.h:705
Utilities::MPI::mean_and_standard_deviation
std::pair< Number, typename numbers::NumberTraits< Number >::real_type > mean_and_standard_deviation(const Iterator begin, const Iterator end, const MPI_Comm &comm)
Utilities::MPI::DuplicatedCommunicator::DuplicatedCommunicator
DuplicatedCommunicator(const MPI_Comm &communicator)
Definition: mpi.h:254
Utilities::MPI::CollectiveMutex::unlock
void unlock(MPI_Comm comm)
Definition: mpi.cc:1138
Utilities::MPI::MinMaxAvg
Definition: mpi.h:687
Utilities::MPI::MPI_InitFinalize::unregister_request
static void unregister_request(MPI_Request &request)
Definition: mpi.cc:911
Utilities::MPI::free_communicator
void free_communicator(MPI_Comm &mpi_communicator)
Definition: mpi.cc:150
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
TransposeTableIterators::Iterator
MatrixTableIterators::Iterator< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Iterator
Definition: table.h:1913
Utilities::MPI::compute_set_union
std::vector< T > compute_set_union(const std::vector< T > &vec, const MPI_Comm &comm)
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
mpi_tags.h
Utilities::MPI::CollectiveMutex::ScopedLock::ScopedLock
ScopedLock(CollectiveMutex &mutex, const MPI_Comm &comm)
Definition: mpi.h:336
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
Utilities::MPI::CollectiveMutex::ScopedLock::comm
const MPI_Comm comm
Definition: mpi.h:359
Tensor
Definition: tensor.h:450
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Utilities::MPI::CollectiveMutex::lock
void lock(MPI_Comm comm)
Definition: mpi.cc:1104
VectorTools::mean
@ mean
Definition: vector_tools_common.h:79
Utilities::MPI::CollectiveMutex
Definition: mpi.h:322
Utilities::MPI::gather
std::vector< T > gather(const MPI_Comm &comm, const T &object_to_send, const unsigned int root_process=0)
Utilities::MPI::create_group
int create_group(const MPI_Comm &comm, const MPI_Group &group, const int tag, MPI_Comm *new_comm)
Definition: mpi.cc:160
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
numbers::NumberTraits::real_type
number real_type
Definition: numbers.h:437
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
Utilities::MPI::CollectiveMutex::~CollectiveMutex
~CollectiveMutex()
Definition: mpi.cc:1091
IndexSet::size_type
types::global_dof_index size_type
Definition: index_set.h:85
Utilities::MPI::DuplicatedCommunicator::comm
MPI_Comm comm
Definition: mpi.h:290
Utilities::MPI::CollectiveMutex::CollectiveMutex
CollectiveMutex()
Definition: mpi.cc:1082
Utilities::MPI::compute_index_owner
std::vector< unsigned int > compute_index_owner(const IndexSet &owned_indices, const IndexSet &indices_to_look_up, const MPI_Comm &comm)
Definition: mpi.cc:1046
unsigned int
Utilities::MPI::MPI_InitFinalize::register_request
static void register_request(MPI_Request &request)
Definition: mpi.cc:902
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Utilities::MPI::this_mpi_process
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
Utilities::MPI::duplicate_communicator
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:139
Utilities::MPI::MinMaxAvg::avg
double avg
Definition: mpi.h:731
Utilities::MPI::MPI_InitFinalize::MPI_InitFinalize
MPI_InitFinalize(int &argc, char **&argv, const unsigned int max_num_threads=numbers::invalid_unsigned_int)
Definition: mpi.cc:741
Utilities::MPI::MinMaxAvg::max_index
unsigned int max_index
Definition: mpi.h:725
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
Utilities::MPI::CollectiveMutex::ScopedLock::~ScopedLock
~ScopedLock()
Definition: mpi.h:346
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Utilities::MPI::CollectiveMutex::locked
bool locked
Definition: mpi.h:394
Utilities::MPI::min_max_avg
MinMaxAvg min_max_avg(const double my_value, const MPI_Comm &mpi_communicator)
Definition: mpi.cc:91
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
Utilities::MPI::compute_point_to_point_communication_pattern
std::vector< unsigned int > compute_point_to_point_communication_pattern(const MPI_Comm &mpi_comm, const std::vector< unsigned int > &destinations)
Definition: mpi.cc:338
Utilities::MPI::MPI_InitFinalize::requests
static std::set< MPI_Request * > requests
Definition: mpi.h:926
config.h
Utilities::MPI::create_ascending_partitioning
std::vector< IndexSet > create_ascending_partitioning(const MPI_Comm &comm, const IndexSet::size_type local_size)
Definition: mpi.cc:242
Utilities::MPI::all_gather
std::vector< T > all_gather(const MPI_Comm &comm, const T &object_to_send)
Utilities::MPI::DuplicatedCommunicator::~DuplicatedCommunicator
~DuplicatedCommunicator()
Definition: mpi.h:266
internal
Definition: aligned_vector.h:369
LAPACKSupport::N
static const char N
Definition: lapack_support.h:159
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Utilities::MPI::MinMaxAvg::min
double min
Definition: mpi.h:699
Utilities::MPI::some_to_some
std::map< unsigned int, T > some_to_some(const MPI_Comm &comm, const std::map< unsigned int, T > &objects_to_send)
Utilities::MPI::job_supports_mpi
bool job_supports_mpi()
Definition: mpi.cc:1030
MPI_Request
numbers.h
Utilities
Definition: cuda.h:31
Utilities::MPI::CollectiveMutex::request
MPI_Request request
Definition: mpi.h:399
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
Utilities::MPI::CollectiveMutex::ScopedLock
Definition: mpi.h:330
int
Utilities::MPI::DuplicatedCommunicator::operator*
const MPI_Comm & operator*() const
Definition: mpi.h:274
Utilities::MPI::MPI_InitFinalize
Definition: mpi.h:828