30 #if !defined(DEAL_II_WITH_MPI) && !defined(DEAL_II_WITH_PETSC)
36 using MPI_Datatype =
int;
38 # ifndef MPI_COMM_WORLD
39 # define MPI_COMM_WORLD 0
41 # ifndef MPI_COMM_SELF
42 # define MPI_COMM_SELF 0
44 # ifndef MPI_REQUEST_NULL
45 # define MPI_REQUEST_NULL 0
73 #ifdef DEAL_II_WITH_MPI
74 # if DEAL_II_MPI_VERSION_GTE(3, 0)
76 # define DEAL_II_MPI_CONST_CAST(expr) (expr)
80 # include <type_traits>
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)
96 template <
int rank,
int dim,
typename Number>
98 template <
int rank,
int dim,
typename Number>
100 template <
typename Number>
121 const unsigned int n_partitions,
176 std::vector<unsigned int>
179 const std::vector<unsigned int> &destinations);
203 const std::vector<unsigned int> &destinations);
431 #ifdef DEAL_II_WITH_MPI
434 const MPI_Group &group,
447 std::vector<IndexSet>
463 #ifdef DEAL_II_WITH_MPI
479 template <
class Iterator,
typename Number =
long double>
480 std::pair<Number, typename numbers::NumberTraits<Number>::real_type>
505 template <
typename T>
518 template <
typename T,
typename U>
520 sum(
const T &values,
const MPI_Comm &mpi_communicator,
U &sums);
531 template <
typename T>
542 template <
int rank,
int dim,
typename Number>
552 template <
int rank,
int dim,
typename Number>
565 template <
typename Number>
590 template <
typename T>
603 template <
typename T,
typename U>
605 max(
const T &values,
const MPI_Comm &mpi_communicator,
U &maxima);
616 template <
typename T>
641 template <
typename T>
654 template <
typename T,
typename U>
656 min(
const T &values,
const MPI_Comm &mpi_communicator,
U &minima);
667 template <
typename T>
762 std::vector<MinMaxAvg>
962 template <
typename T>
963 std::map<unsigned int, T>
965 const std::map<unsigned int, T> &objects_to_send);
982 template <
typename T>
1003 template <
typename T>
1006 const T & object_to_send,
1007 const unsigned int root_process = 0);
1053 std::vector<unsigned int>
1055 const IndexSet &indices_to_look_up,
1065 template <
typename T>
1072 template <
typename T>
1080 template <
typename T>
1082 all_reduce(
const MPI_Op & mpi_op,
1089 template <
typename T,
unsigned int N>
1091 sum(
const T (&values)[
N],
const MPI_Comm &mpi_communicator,
T (&sums)[
N])
1093 internal::all_reduce(MPI_SUM,
1099 template <
typename T,
unsigned int N>
1101 max(
const T (&values)[
N],
const MPI_Comm &mpi_communicator,
T (&maxima)[
N])
1103 internal::all_reduce(MPI_MAX,
1109 template <
typename T,
unsigned int N>
1111 min(
const T (&values)[
N],
const MPI_Comm &mpi_communicator,
T (&minima)[
N])
1113 internal::all_reduce(MPI_MIN,
1119 template <
typename T>
1120 std::map<unsigned int, T>
1122 const std::map<unsigned int, T> &objects_to_send)
1124 # ifndef DEAL_II_WITH_MPI
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;
1135 std::map<unsigned int, T> received_objects;
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;
1143 send_to.emplace_back(m.first);
1145 const unsigned int n_point_point_communications =
1149 static CollectiveMutex mutex;
1150 CollectiveMutex::ScopedLock lock(mutex, comm);
1156 if (send_to.size() == 0 && n_point_point_communications == 0)
1157 return received_objects;
1163 std::vector<std::vector<char>> buffers_to_send(send_to.size());
1164 std::vector<MPI_Request> buffer_send_requests(send_to.size());
1167 for (
const auto &rank_obj : objects_to_send)
1168 if (rank_obj.first != my_proc)
1170 const auto &rank = rank_obj.first;
1172 const int ierr = MPI_Isend(buffers_to_send[i].data(),
1173 buffers_to_send[i].size(),
1178 &buffer_send_requests[i]);
1186 std::vector<char> buffer;
1188 for (
unsigned int i = 0; i < n_point_point_communications; ++i)
1192 int ierr = MPI_Probe(MPI_ANY_SOURCE, mpi_tag, comm, &status);
1197 ierr = MPI_Get_count(&status, MPI_CHAR, &len);
1202 const unsigned int rank = status.MPI_SOURCE;
1205 ierr = MPI_Recv(buffer.data(),
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);
1221 const int ierr = MPI_Waitall(send_to.size(),
1222 buffer_send_requests.data(),
1223 MPI_STATUSES_IGNORE);
1226 return received_objects;
1227 # endif // deal.II with MPI
1230 template <
typename T>
1237 # ifndef DEAL_II_WITH_MPI
1239 std::vector<T> v(1,
object);
1246 int n_local_data = buffer.size();
1249 std::vector<int> size_all_data(n_procs, 0);
1253 &n_local_data, 1, MPI_INT, size_all_data.data(), 1, MPI_INT, comm);
1257 std::vector<int> rdispls(n_procs);
1259 for (
unsigned int i = 1; i < n_procs; ++i)
1260 rdispls[i] = rdispls[i - 1] + size_all_data[i - 1];
1263 std::vector<char> received_unrolled_buffer(rdispls.back() +
1264 size_all_data.back());
1266 MPI_Allgatherv(buffer.data(),
1269 received_unrolled_buffer.data(),
1270 size_all_data.data(),
1275 std::vector<T> received_objects(n_procs);
1276 for (
unsigned int i = 0; i < n_procs; ++i)
1278 std::vector<char> local_buffer(received_unrolled_buffer.begin() +
1280 received_unrolled_buffer.begin() +
1281 rdispls[i] + size_all_data[i]);
1282 received_objects[i] = Utilities::unpack<T>(local_buffer);
1285 return received_objects;
1289 template <
typename T>
1292 const T & object_to_send,
1293 const unsigned int root_process)
1295 # ifndef DEAL_II_WITH_MPI
1298 std::vector<T> v(1, object_to_send);
1307 int n_local_data = buffer.size();
1311 std::vector<int> size_all_data;
1312 if (my_rank == root_process)
1313 size_all_data.resize(n_procs, 0);
1316 int ierr = MPI_Gather(&n_local_data,
1319 size_all_data.data(),
1328 std::vector<int> rdispls;
1329 if (my_rank == root_process)
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];
1336 std::vector<char> received_unrolled_buffer;
1337 if (my_rank == root_process)
1338 received_unrolled_buffer.resize(rdispls.back() + size_all_data.back());
1340 ierr = MPI_Gatherv(buffer.data(),
1343 received_unrolled_buffer.data(),
1344 size_all_data.data(),
1351 std::vector<T> received_objects;
1353 if (my_rank == root_process)
1355 received_objects.resize(n_procs);
1357 for (
unsigned int i = 0; i < n_procs; ++i)
1359 const std::vector<char> local_buffer(
1360 received_unrolled_buffer.begin() + rdispls[i],
1361 received_unrolled_buffer.begin() + rdispls[i] +
1363 received_objects[i] = Utilities::unpack<T>(local_buffer);
1366 return received_objects;
1371 # ifdef DEAL_II_WITH_MPI
1372 template <
class Iterator,
typename Number>
1373 std::pair<Number, typename numbers::NumberTraits<Number>::real_type>
1384 const Number
sum = std::accumulate(
begin,
end, Number(0.));
1391 std::for_each(
begin,
end, [&
mean, &sq_sum](
const Number &v) {
1395 return std::make_pair(
mean,
1396 std::sqrt(sq_sum /
static_cast<Std
>(size - 1)));