14#include <dolfinx/common/Timer.h>
15#include <dolfinx/common/log.h>
16#include <dolfinx/graph/AdjacencyList.h>
25#define MPICH_IGNORE_CXX_SEEK 1
45 explicit Comm(MPI_Comm
comm,
bool duplicate =
true);
63 MPI_Comm
comm()
const noexcept;
71int rank(MPI_Comm comm);
75int size(MPI_Comm comm);
98 const std::int64_t n = N /
size;
99 const std::int64_t r = N %
size;
103 return {
rank * (n + 1),
rank * (n + 1) + n + 1};
105 return {
rank * n + r,
rank * n + r + n};
119 const std::size_t n = N /
size;
120 const std::size_t r = N %
size;
122 if (index < r * (n + 1))
125 return index / (n + 1);
130 return r + (index - r * (n + 1)) / n;
158 std::span<const int> edges);
185 std::span<const int> edges);
207std::pair<std::vector<std::int32_t>, std::vector<T>>
209 std::array<std::int64_t, 2> shape,
210 std::int64_t rank_offset);
235 std::span<const std::int64_t> indices,
236 std::span<const T> x,
237 std::array<std::int64_t, 2> shape,
238 std::int64_t rank_offset);
265 std::span<const std::int64_t> indices,
266 std::span<const T> x,
int shape1);
277 if constexpr (std::is_same_v<T, float>)
279 else if constexpr (std::is_same_v<T, double>)
281 else if constexpr (std::is_same_v<T, std::complex<double>>)
282 return MPI_C_DOUBLE_COMPLEX;
283 else if constexpr (std::is_same_v<T, std::complex<float>>)
284 return MPI_C_FLOAT_COMPLEX;
285 else if constexpr (std::is_same_v<T, short int>)
287 else if constexpr (std::is_same_v<T, int>)
289 else if constexpr (std::is_same_v<T, unsigned int>)
291 else if constexpr (std::is_same_v<T, long int>)
293 else if constexpr (std::is_same_v<T, unsigned long>)
294 return MPI_UNSIGNED_LONG;
295 else if constexpr (std::is_same_v<T, long long>)
296 return MPI_LONG_LONG;
297 else if constexpr (std::is_same_v<T, unsigned long long>)
298 return MPI_UNSIGNED_LONG_LONG;
299 else if constexpr (std::is_same_v<T, bool>)
301 else if constexpr (std::is_same_v<T, std::int8_t>)
305 static_assert(!std::is_same_v<T, T>);
310std::pair<std::vector<std::int32_t>, std::vector<T>>
312 std::array<std::int64_t, 2> shape,
313 std::int64_t rank_offset)
317 assert(x.size() % shape[1] == 0);
318 const std::int32_t shape0_local = x.size() / shape[1];
320 LOG(2) <<
"Sending data to post offices (distribute_to_postoffice)";
323 std::vector<int> row_to_dest(shape0_local);
324 for (std::int32_t i = 0; i < shape0_local; ++i)
327 row_to_dest[i] = dest;
332 std::vector<std::array<std::int32_t, 2>> dest_to_index;
333 dest_to_index.reserve(shape0_local);
334 for (std::int32_t i = 0; i < shape0_local; ++i)
336 std::size_t idx = i + rank_offset;
338 dest_to_index.push_back({dest, i});
340 std::sort(dest_to_index.begin(), dest_to_index.end());
344 std::vector<int> dest;
345 std::vector<std::int32_t> num_items_per_dest,
346 pos_to_neigh_rank(shape0_local, -1);
348 auto it = dest_to_index.begin();
349 while (it != dest_to_index.end())
351 const int neigh_rank = dest.size();
354 dest.push_back((*it)[0]);
358 = std::find_if(it, dest_to_index.end(),
359 [r = dest.back()](
auto& idx) { return idx[0] != r; });
362 num_items_per_dest.push_back(std::distance(it, it1));
365 for (
auto e = it; e != it1; ++e)
366 pos_to_neigh_rank[(*e)[1]] = neigh_rank;
376 <<
"Number of neighbourhood source ranks in distribute_to_postoffice: "
381 int err = MPI_Dist_graph_create_adjacent(
382 comm, src.size(), src.data(), MPI_UNWEIGHTED, dest.size(), dest.data(),
383 MPI_UNWEIGHTED, MPI_INFO_NULL,
false, &neigh_comm);
387 std::vector<std::int32_t> send_disp = {0};
388 std::partial_sum(num_items_per_dest.begin(), num_items_per_dest.end(),
389 std::back_inserter(send_disp));
392 std::vector<T> send_buffer_data(shape[1] * send_disp.back());
393 std::vector<std::int64_t> send_buffer_index(send_disp.back());
395 std::vector<std::int32_t> send_offsets = send_disp;
396 for (std::int32_t i = 0; i < shape0_local; ++i)
398 if (
int neigh_dest = pos_to_neigh_rank[i]; neigh_dest != -1)
400 std::size_t pos = send_offsets[neigh_dest];
401 send_buffer_index[pos] = i + rank_offset;
402 std::copy_n(std::next(x.begin(), i * shape[1]), shape[1],
403 std::next(send_buffer_data.begin(), shape[1] * pos));
404 ++send_offsets[neigh_dest];
411 std::vector<int> num_items_recv(src.size());
412 num_items_per_dest.reserve(1);
413 num_items_recv.reserve(1);
414 err = MPI_Neighbor_alltoall(num_items_per_dest.data(), 1, MPI_INT,
415 num_items_recv.data(), 1, MPI_INT, neigh_comm);
419 std::vector<std::int32_t> recv_disp(num_items_recv.size() + 1, 0);
420 std::partial_sum(num_items_recv.begin(), num_items_recv.end(),
421 std::next(recv_disp.begin()));
424 std::vector<std::int64_t> recv_buffer_index(recv_disp.back());
425 err = MPI_Neighbor_alltoallv(
426 send_buffer_index.data(), num_items_per_dest.data(), send_disp.data(),
427 MPI_INT64_T, recv_buffer_index.data(), num_items_recv.data(),
428 recv_disp.data(), MPI_INT64_T, neigh_comm);
432 MPI_Datatype compound_type;
433 MPI_Type_contiguous(shape[1], dolfinx::MPI::mpi_type<T>(), &compound_type);
434 MPI_Type_commit(&compound_type);
435 std::vector<T> recv_buffer_data(shape[1] * recv_disp.back());
436 err = MPI_Neighbor_alltoallv(
437 send_buffer_data.data(), num_items_per_dest.data(), send_disp.data(),
438 compound_type, recv_buffer_data.data(), num_items_recv.data(),
439 recv_disp.data(), compound_type, neigh_comm);
441 err = MPI_Type_free(&compound_type);
443 err = MPI_Comm_free(&neigh_comm);
446 LOG(2) <<
"Completed send data to post offices.";
450 std::vector<std::int32_t> index_local(recv_buffer_index.size());
451 std::transform(recv_buffer_index.cbegin(), recv_buffer_index.cend(),
452 index_local.begin(), [r0](
auto idx) { return idx - r0; });
454 return {index_local, recv_buffer_data};
459 std::span<const std::int64_t> indices,
460 std::span<const T> x,
461 std::array<std::int64_t, 2> shape,
462 std::int64_t rank_offset)
465 assert(shape[1] > 0);
469 assert(x.size() % shape[1] == 0);
470 const std::int64_t shape0_local = x.size() / shape[1];
477 comm, x, {shape[0], shape[1]}, rank_offset);
478 assert(post_indices.size() == post_x.size() / shape[1]);
484 std::vector<std::tuple<int, std::int64_t, std::int32_t>> src_to_index;
485 for (std::size_t i = 0; i < indices.size(); ++i)
487 std::size_t idx = indices[i];
489 src_to_index.push_back({src, idx, i});
491 std::sort(src_to_index.begin(), src_to_index.end());
495 std::vector<std::int32_t> num_items_per_src;
496 std::vector<int> src;
498 auto it = src_to_index.begin();
499 while (it != src_to_index.end())
501 src.push_back(std::get<0>(*it));
502 auto it1 = std::find_if(it, src_to_index.end(),
503 [r = src.back()](
auto& idx)
504 { return std::get<0>(idx) != r; });
505 num_items_per_src.push_back(std::distance(it, it1));
512 const std::vector<int> dest
514 LOG(INFO) <<
"Neighbourhood destination ranks from post office in "
515 "distribute_data (rank, num dests, num dests/mpi_size): "
516 <<
rank <<
", " << dest.size() <<
", "
517 <<
static_cast<double>(dest.size()) /
size;
521 MPI_Comm neigh_comm0;
522 int err = MPI_Dist_graph_create_adjacent(
523 comm, dest.size(), dest.data(), MPI_UNWEIGHTED, src.size(), src.data(),
524 MPI_UNWEIGHTED, MPI_INFO_NULL,
false, &neigh_comm0);
528 std::vector<int> num_items_recv(dest.size());
529 num_items_per_src.reserve(1);
530 num_items_recv.reserve(1);
531 err = MPI_Neighbor_alltoall(num_items_per_src.data(), 1, MPI_INT,
532 num_items_recv.data(), 1, MPI_INT, neigh_comm0);
536 std::vector<std::int32_t> send_disp = {0};
537 std::partial_sum(num_items_per_src.begin(), num_items_per_src.end(),
538 std::back_inserter(send_disp));
539 std::vector<std::int32_t> recv_disp = {0};
540 std::partial_sum(num_items_recv.begin(), num_items_recv.end(),
541 std::back_inserter(recv_disp));
545 assert(send_disp.back() == (
int)src_to_index.size());
546 std::vector<std::int64_t> send_buffer_index(src_to_index.size());
547 std::transform(src_to_index.cbegin(), src_to_index.cend(),
548 send_buffer_index.begin(),
549 [](
auto& x) { return std::get<1>(x); });
552 std::vector<std::int64_t> recv_buffer_index(recv_disp.back());
553 err = MPI_Neighbor_alltoallv(
554 send_buffer_index.data(), num_items_per_src.data(), send_disp.data(),
555 MPI_INT64_T, recv_buffer_index.data(), num_items_recv.data(),
556 recv_disp.data(), MPI_INT64_T, neigh_comm0);
559 err = MPI_Comm_free(&neigh_comm0);
568 const std::array<std::int64_t, 2> postoffice_range
570 std::vector<std::int32_t> post_indices_map(
571 postoffice_range[1] - postoffice_range[0], -1);
572 for (std::size_t i = 0; i < post_indices.size(); ++i)
574 assert(post_indices[i] < (
int)post_indices_map.size());
575 post_indices_map[post_indices[i]] = i;
579 std::vector<T> send_buffer_data(shape[1] * recv_disp.back());
580 for (std::size_t p = 0; p < recv_disp.size() - 1; ++p)
582 int offset = recv_disp[p];
583 for (std::int32_t i = recv_disp[p]; i < recv_disp[p + 1]; ++i)
585 std::int64_t index = recv_buffer_index[i];
586 if (index >= rank_offset and index < (rank_offset + shape0_local))
589 std::int32_t local_index = index - rank_offset;
590 std::copy_n(std::next(x.begin(), shape[1] * local_index), shape[1],
591 std::next(send_buffer_data.begin(), shape[1] * offset));
596 auto local_index = index - postoffice_range[0];
597 std::int32_t pos = post_indices_map[local_index];
599 std::copy_n(std::next(post_x.begin(), shape[1] * pos), shape[1],
600 std::next(send_buffer_data.begin(), shape[1] * offset));
607 err = MPI_Dist_graph_create_adjacent(
608 comm, src.size(), src.data(), MPI_UNWEIGHTED, dest.size(), dest.data(),
609 MPI_UNWEIGHTED, MPI_INFO_NULL,
false, &neigh_comm0);
612 MPI_Datatype compound_type0;
613 MPI_Type_contiguous(shape[1], dolfinx::MPI::mpi_type<T>(), &compound_type0);
614 MPI_Type_commit(&compound_type0);
616 std::vector<T> recv_buffer_data(shape[1] * send_disp.back());
617 err = MPI_Neighbor_alltoallv(
618 send_buffer_data.data(), num_items_recv.data(), recv_disp.data(),
619 compound_type0, recv_buffer_data.data(), num_items_per_src.data(),
620 send_disp.data(), compound_type0, neigh_comm0);
623 err = MPI_Type_free(&compound_type0);
625 err = MPI_Comm_free(&neigh_comm0);
628 std::vector<std::int32_t> index_pos_to_buffer(indices.size(), -1);
629 for (std::size_t i = 0; i < src_to_index.size(); ++i)
630 index_pos_to_buffer[std::get<2>(src_to_index[i])] = i;
633 std::vector<T> x_new(shape[1] * indices.size());
634 for (std::size_t i = 0; i < indices.size(); ++i)
636 const std::int64_t index = indices[i];
637 if (index >= rank_offset and index < (rank_offset + shape0_local))
640 auto local_index = index - rank_offset;
641 std::copy_n(std::next(x.begin(), shape[1] * local_index), shape[1],
642 std::next(x_new.begin(), shape[1] * i));
649 auto local_index = index - postoffice_range[0];
650 std::int32_t pos = post_indices_map[local_index];
652 std::copy_n(std::next(post_x.begin(), shape[1] * pos), shape[1],
653 std::next(x_new.begin(), shape[1] * i));
658 std::int32_t pos = index_pos_to_buffer[i];
660 std::copy_n(std::next(recv_buffer_data.begin(), shape[1] * pos),
661 shape[1], std::next(x_new.begin(), shape[1] * i));
671 std::span<const std::int64_t> indices,
672 std::span<const T> x,
int shape1)
675 assert(x.size() % shape1 == 0);
676 const std::int64_t shape0_local = x.size() / shape1;
678 std::int64_t shape0(0), rank_offset(0);
680 = MPI_Allreduce(&shape0_local, &shape0, 1, MPI_INT64_T, MPI_SUM, comm);
682 err = MPI_Exscan(&shape0_local, &rank_offset, 1, MPI_INT64_T, MPI_SUM, comm);
A duplicate MPI communicator and manage lifetime of the communicator.
Definition: MPI.h:42
~Comm()
Destructor (frees wrapped communicator)
Definition: MPI.cpp:35
MPI_Comm comm() const noexcept
Return the underlying MPI_Comm object.
Definition: MPI.cpp:61
A timer can be used for timing tasks. The basic usage is.
Definition: Timer.h:31
MPI support functionality.
Definition: MPI.h:30
std::pair< std::vector< std::int32_t >, std::vector< T > > distribute_to_postoffice(MPI_Comm comm, std::span< const T > x, std::array< std::int64_t, 2 > shape, std::int64_t rank_offset)
Distribute row data to 'post office' ranks.
Definition: MPI.h:311
std::vector< int > compute_graph_edges_nbx(MPI_Comm comm, std::span< const int > edges)
Determine incoming graph edges using the NBX consensus algorithm.
Definition: MPI.cpp:163
std::vector< T > distribute_from_postoffice(MPI_Comm comm, std::span< const std::int64_t > indices, std::span< const T > x, std::array< std::int64_t, 2 > shape, std::int64_t rank_offset)
Distribute rows of a rectangular data array from post office ranks to ranks where they are required.
Definition: MPI.h:458
constexpr int index_owner(int size, std::size_t index, std::size_t N)
Return which rank owns index in global range [0, N - 1] (inverse of MPI::local_range).
Definition: MPI.h:114
std::vector< T > distribute_data(MPI_Comm comm, std::span< const std::int64_t > indices, std::span< const T > x, int shape1)
Distribute rows of a rectangular data array to ranks where they are required (scalable version).
Definition: MPI.h:670
std::vector< int > compute_graph_edges_pcx(MPI_Comm comm, std::span< const int > edges)
Determine incoming graph edges using the PCX consensus algorithm.
Definition: MPI.cpp:96
void check_error(MPI_Comm comm, int code)
Check MPI error code. If the error code is not equal to MPI_SUCCESS, then std::abort is called.
Definition: MPI.cpp:79
int size(MPI_Comm comm)
Return size of the group (number of processes) associated with the communicator.
Definition: MPI.cpp:71
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition: MPI.cpp:63
constexpr std::array< std::int64_t, 2 > local_range(int rank, std::int64_t N, int size)
Return local range for the calling process, partitioning the global [0, N - 1] range across all ranks...
Definition: MPI.h:90
constexpr MPI_Datatype mpi_type()
MPI Type.
Definition: MPI.h:275
tag
MPI communication tags.
Definition: MPI.h:34