DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
MPI.h
1// Copyright (C) 2007-2023 Magnus Vikstrøm, Garth N. Wells and Paul T. Kühner
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "Timer.h"
10#include "log.h"
11#include "types.h"
12#include <algorithm>
13#include <array>
14#include <cassert>
15#include <complex>
16#include <cstdint>
17#include <numeric>
18#include <set>
19#include <span>
20#include <tuple>
21#include <type_traits>
22#include <utility>
23#include <vector>
24
25#define MPICH_IGNORE_CXX_SEEK 1
26#include <mpi.h>
27
29namespace dolfinx::MPI
30{
32enum class tag : std::uint16_t
33{
34 consensus_pcx = 1200,
35 consensus_pex = 1201,
36 consensus_nbx = 1202,
37};
38
41class Comm
42{
43public:
45 explicit Comm(MPI_Comm comm, bool duplicate = true);
46
48 Comm(const Comm& comm) noexcept;
49
51 Comm(Comm&& comm) noexcept;
52
53 // Disable copy assignment operator
54 Comm& operator=(const Comm& comm) = delete;
55
57 Comm& operator=(Comm&& comm) noexcept;
58
60 ~Comm();
61
63 MPI_Comm comm() const noexcept;
64
65private:
66 // MPI communicator
67 MPI_Comm _comm;
68};
69
71int rank(MPI_Comm comm);
72
75int size(MPI_Comm comm);
76
81void check_error(MPI_Comm comm, int code);
82
89constexpr std::array<std::int64_t, 2> local_range(int rank, std::int64_t N,
90 int size)
91{
92 assert(rank >= 0);
93 assert(N >= 0);
94 assert(size > 0);
95
96 // Compute number of items per rank and remainder
97 const std::int64_t n = N / size;
98 const std::int64_t r = N % size;
99
100 // Compute local range
101 if (rank < r)
102 return {rank * (n + 1), rank * (n + 1) + n + 1};
103 else
104 return {rank * n + r, rank * n + r + n};
105}
106
113constexpr int index_owner(int size, std::size_t index, std::size_t N)
114{
115 assert(index < N);
116
117 // Compute number of items per rank and remainder
118 const std::size_t n = N / size;
119 const std::size_t r = N % size;
120
121 if (index < r * (n + 1))
122 {
123 // First r ranks own n + 1 indices
124 return index / (n + 1);
125 }
126 else
127 {
128 // Remaining ranks own n indices
129 return r + (index - r * (n + 1)) / n;
130 }
131}
132
156std::vector<int> compute_graph_edges_pcx(MPI_Comm comm,
157 std::span<const int> edges);
158
184std::vector<int>
185compute_graph_edges_nbx(MPI_Comm comm, std::span<const int> edges,
186 int tag = static_cast<int>(tag::consensus_nbx));
187
207template <typename U>
208std::pair<std::vector<std::int32_t>,
209 std::vector<typename std::remove_reference_t<typename U::value_type>>>
210distribute_to_postoffice(MPI_Comm comm, const U& x,
211 std::array<std::int64_t, 2> shape,
212 std::int64_t rank_offset);
213
235template <typename U>
236std::vector<typename std::remove_reference_t<typename U::value_type>>
237distribute_from_postoffice(MPI_Comm comm, std::span<const std::int64_t> indices,
238 const U& x, std::array<std::int64_t, 2> shape,
239 std::int64_t rank_offset);
240
261template <typename U>
262std::vector<typename std::remove_reference_t<typename U::value_type>>
263distribute_data(MPI_Comm comm0, std::span<const std::int64_t> indices,
264 MPI_Comm comm1, const U& x, int shape1);
265
266template <typename T>
267struct dependent_false : std::false_type
268{
269};
270
272
274template <typename T>
276
279template <typename T>
281
284#define MAP_TO_MPI_TYPE(cpp_t, mpi_t) \
285 template <> \
286 struct mpi_type_mapping<cpp_t> \
287 { \
288 static inline MPI_Datatype type = mpi_t; \
289 };
290
294MAP_TO_MPI_TYPE(float, MPI_FLOAT)
295MAP_TO_MPI_TYPE(double, MPI_DOUBLE)
296MAP_TO_MPI_TYPE(std::complex<float>, MPI_C_FLOAT_COMPLEX)
297MAP_TO_MPI_TYPE(std::complex<double>, MPI_C_DOUBLE_COMPLEX)
298MAP_TO_MPI_TYPE(std::int8_t, MPI_INT8_T)
299MAP_TO_MPI_TYPE(std::int16_t, MPI_INT16_T)
300MAP_TO_MPI_TYPE(std::int32_t, MPI_INT32_T)
301MAP_TO_MPI_TYPE(std::int64_t, MPI_INT64_T)
302MAP_TO_MPI_TYPE(std::uint8_t, MPI_UINT8_T)
303MAP_TO_MPI_TYPE(std::uint16_t, MPI_UINT16_T)
304MAP_TO_MPI_TYPE(std::uint32_t, MPI_UINT32_T)
305MAP_TO_MPI_TYPE(std::uint64_t, MPI_UINT64_T)
308
309//---------------------------------------------------------------------------
310template <typename U>
311std::pair<std::vector<std::int32_t>,
312 std::vector<typename std::remove_reference_t<typename U::value_type>>>
313distribute_to_postoffice(MPI_Comm comm, const U& x,
314 std::array<std::int64_t, 2> shape,
315 std::int64_t rank_offset)
316{
317 assert(rank_offset >= 0 or x.empty());
318 using T = typename std::remove_reference_t<typename U::value_type>;
319
320 const int size = dolfinx::MPI::size(comm);
321 const int rank = dolfinx::MPI::rank(comm);
322 assert(x.size() % shape[1] == 0);
323 const std::int32_t shape0_local = x.size() / shape[1];
324
325 spdlog::debug("Sending data to post offices (distribute_to_postoffice)");
326
327 // Post office ranks will receive data from this rank
328 std::vector<int> row_to_dest(shape0_local);
329 for (std::int32_t i = 0; i < shape0_local; ++i)
330 {
331 int dest = MPI::index_owner(size, i + rank_offset, shape[0]);
332 row_to_dest[i] = dest;
333 }
334
335 // Build list of (dest, positions) for each row that doesn't belong to
336 // this rank, then sort
337 std::vector<std::array<std::int32_t, 2>> dest_to_index;
338 dest_to_index.reserve(shape0_local);
339 for (std::int32_t i = 0; i < shape0_local; ++i)
340 {
341 std::size_t idx = i + rank_offset;
342 if (int dest = MPI::index_owner(size, idx, shape[0]); dest != rank)
343 dest_to_index.push_back({dest, i});
344 }
345 std::ranges::sort(dest_to_index);
346
347 // Build list of neighbour src ranks and count number of items (rows
348 // of x) to receive from each src post office (by neighbourhood rank)
349 std::vector<int> dest;
350 std::vector<std::int32_t> num_items_per_dest,
351 pos_to_neigh_rank(shape0_local, -1);
352 {
353 auto it = dest_to_index.begin();
354 while (it != dest_to_index.end())
355 {
356 const int neigh_rank = dest.size();
357
358 // Store global rank
359 dest.push_back((*it)[0]);
360
361 // Find iterator to next global rank
362 auto it1
363 = std::find_if(it, dest_to_index.end(),
364 [r = dest.back()](auto& idx) { return idx[0] != r; });
365
366 // Store number of items for current rank
367 num_items_per_dest.push_back(std::distance(it, it1));
368
369 // Map from local x index to local destination rank
370 for (auto e = it; e != it1; ++e)
371 pos_to_neigh_rank[(*e)[1]] = neigh_rank;
372
373 // Advance iterator
374 it = it1;
375 }
376 }
377
378 // Determine source ranks
379 const std::vector<int> src = MPI::compute_graph_edges_nbx(comm, dest);
380 spdlog::info(
381 "Number of neighbourhood source ranks in distribute_to_postoffice: {}",
382 static_cast<int>(src.size()));
383
384 // Create neighbourhood communicator for sending data to post offices
385 MPI_Comm neigh_comm;
386 int err = MPI_Dist_graph_create_adjacent(
387 comm, src.size(), src.data(), MPI_UNWEIGHTED, dest.size(), dest.data(),
388 MPI_UNWEIGHTED, MPI_INFO_NULL, false, &neigh_comm);
389 dolfinx::MPI::check_error(comm, err);
390
391 // Compute send displacements
392 std::vector<std::int32_t> send_disp{0};
393 std::partial_sum(num_items_per_dest.begin(), num_items_per_dest.end(),
394 std::back_inserter(send_disp));
395
396 // Pack send buffers
397 std::vector<T> send_buffer_data(shape[1] * send_disp.back());
398 std::vector<std::int64_t> send_buffer_index(send_disp.back());
399 {
400 std::vector<std::int32_t> send_offsets = send_disp;
401 for (std::int32_t i = 0; i < shape0_local; ++i)
402 {
403 if (int neigh_dest = pos_to_neigh_rank[i]; neigh_dest != -1)
404 {
405 std::size_t pos = send_offsets[neigh_dest];
406 send_buffer_index[pos] = i + rank_offset;
407 std::copy_n(std::next(x.begin(), i * shape[1]), shape[1],
408 std::next(send_buffer_data.begin(), shape[1] * pos));
409 ++send_offsets[neigh_dest];
410 }
411 }
412 }
413
414 // Send number of items to post offices (destination) that I will be
415 // sending
416 std::vector<int> num_items_recv(src.size());
417 num_items_per_dest.reserve(1);
418 num_items_recv.reserve(1);
419 err = MPI_Neighbor_alltoall(num_items_per_dest.data(), 1, MPI_INT,
420 num_items_recv.data(), 1, MPI_INT, neigh_comm);
421 dolfinx::MPI::check_error(comm, err);
422
423 // Prepare receive displacement and buffers
424 std::vector<std::int32_t> recv_disp(num_items_recv.size() + 1, 0);
425 std::partial_sum(num_items_recv.begin(), num_items_recv.end(),
426 std::next(recv_disp.begin()));
427
428 // Send/receive global indices
429 std::vector<std::int64_t> recv_buffer_index(recv_disp.back());
430 err = MPI_Neighbor_alltoallv(
431 send_buffer_index.data(), num_items_per_dest.data(), send_disp.data(),
432 MPI_INT64_T, recv_buffer_index.data(), num_items_recv.data(),
433 recv_disp.data(), MPI_INT64_T, neigh_comm);
434 dolfinx::MPI::check_error(comm, err);
435
436 // Send/receive data (x)
437 MPI_Datatype compound_type;
438 MPI_Type_contiguous(shape[1], dolfinx::MPI::mpi_t<T>, &compound_type);
439 MPI_Type_commit(&compound_type);
440 std::vector<T> recv_buffer_data(shape[1] * recv_disp.back());
441 err = MPI_Neighbor_alltoallv(
442 send_buffer_data.data(), num_items_per_dest.data(), send_disp.data(),
443 compound_type, recv_buffer_data.data(), num_items_recv.data(),
444 recv_disp.data(), compound_type, neigh_comm);
445 dolfinx::MPI::check_error(comm, err);
446 err = MPI_Type_free(&compound_type);
447 dolfinx::MPI::check_error(comm, err);
448 err = MPI_Comm_free(&neigh_comm);
449 dolfinx::MPI::check_error(comm, err);
450
451 spdlog::debug("Completed send data to post offices.");
452
453 // Convert to local indices
454 const std::int64_t r0 = MPI::local_range(rank, shape[0], size)[0];
455 std::vector<std::int32_t> index_local(recv_buffer_index.size());
456 std::ranges::transform(recv_buffer_index, index_local.begin(),
457 [r0](auto idx) { return idx - r0; });
458
459 return {index_local, recv_buffer_data};
460}
461//---------------------------------------------------------------------------
462template <typename U>
463std::vector<typename std::remove_reference_t<typename U::value_type>>
464distribute_from_postoffice(MPI_Comm comm, std::span<const std::int64_t> indices,
465 const U& x, std::array<std::int64_t, 2> shape,
466 std::int64_t rank_offset)
467{
468 assert(rank_offset >= 0 or x.empty());
469 using T = typename std::remove_reference_t<typename U::value_type>;
470
471 common::Timer timer("Distribute row-wise data (scalable)");
472 assert(shape[1] > 0);
473
474 const int size = dolfinx::MPI::size(comm);
475 const int rank = dolfinx::MPI::rank(comm);
476 assert(x.size() % shape[1] == 0);
477 const std::int64_t shape0_local = x.size() / shape[1];
478
479 // 0. Send x data to/from post offices
480
481 // Send receive x data to post office (only for rows that need to be
482 // communicated)
483 auto [post_indices, post_x] = dolfinx::MPI::distribute_to_postoffice(
484 comm, x, {shape[0], shape[1]}, rank_offset);
485 assert(post_indices.size() == post_x.size() / shape[1]);
486
487 // 1. Send request to post office ranks for data
488
489 // Build list of (src, global index, global, index position) for each
490 // entry in 'indices' that doesn't belong to this rank, then sort
491 std::vector<std::tuple<int, std::int64_t, std::int32_t>> src_to_index;
492 for (std::size_t i = 0; i < indices.size(); ++i)
493 {
494 std::size_t idx = indices[i];
495 if (int src = dolfinx::MPI::index_owner(size, idx, shape[0]); src != rank)
496 src_to_index.push_back({src, idx, i});
497 }
498 std::ranges::sort(src_to_index);
499
500 // Build list is neighbour src ranks and count number of items (rows
501 // of x) to receive from each src post office (by neighbourhood rank)
502 std::vector<std::int32_t> num_items_per_src;
503 std::vector<int> src;
504 {
505 auto it = src_to_index.begin();
506 while (it != src_to_index.end())
507 {
508 src.push_back(std::get<0>(*it));
509 auto it1
510 = std::find_if(it, src_to_index.end(), [r = src.back()](auto& idx)
511 { return std::get<0>(idx) != r; });
512 num_items_per_src.push_back(std::distance(it, it1));
513 it = it1;
514 }
515 }
516
517 // Determine 'delivery' destination ranks (ranks that want data from
518 // me)
519 const std::vector<int> dest
521 spdlog::info(
522 "Neighbourhood destination ranks from post office in "
523 "distribute_data (rank, num dests, num dests/mpi_size): {}, {}, {}",
524 rank, static_cast<int>(dest.size()),
525 static_cast<double>(dest.size()) / size);
526
527 // Create neighbourhood communicator for sending data to post offices
528 // (src), and receiving data form my send my post office
529 MPI_Comm neigh_comm0;
530 int err = MPI_Dist_graph_create_adjacent(
531 comm, dest.size(), dest.data(), MPI_UNWEIGHTED, src.size(), src.data(),
532 MPI_UNWEIGHTED, MPI_INFO_NULL, false, &neigh_comm0);
533 dolfinx::MPI::check_error(comm, err);
534
535 // Communicate number of requests to each source
536 std::vector<int> num_items_recv(dest.size());
537 num_items_per_src.reserve(1);
538 num_items_recv.reserve(1);
539 err = MPI_Neighbor_alltoall(num_items_per_src.data(), 1, MPI_INT,
540 num_items_recv.data(), 1, MPI_INT, neigh_comm0);
541 dolfinx::MPI::check_error(comm, err);
542
543 // Prepare send/receive displacements
544 std::vector<std::int32_t> send_disp{0};
545 std::partial_sum(num_items_per_src.begin(), num_items_per_src.end(),
546 std::back_inserter(send_disp));
547 std::vector<std::int32_t> recv_disp = {0};
548 std::partial_sum(num_items_recv.begin(), num_items_recv.end(),
549 std::back_inserter(recv_disp));
550
551 // Pack my requested indices (global) in send buffer ready to send to
552 // post offices
553 assert(send_disp.back() == (int)src_to_index.size());
554 std::vector<std::int64_t> send_buffer_index(src_to_index.size());
555 std::ranges::transform(src_to_index, send_buffer_index.begin(),
556 [](auto x) { return std::get<1>(x); });
557
558 // Prepare the receive buffer
559 std::vector<std::int64_t> recv_buffer_index(recv_disp.back());
560 err = MPI_Neighbor_alltoallv(
561 send_buffer_index.data(), num_items_per_src.data(), send_disp.data(),
562 MPI_INT64_T, recv_buffer_index.data(), num_items_recv.data(),
563 recv_disp.data(), MPI_INT64_T, neigh_comm0);
564 dolfinx::MPI::check_error(comm, err);
565
566 err = MPI_Comm_free(&neigh_comm0);
567 dolfinx::MPI::check_error(comm, err);
568
569 // 2. Send data (rows of x) from post office back to requesting ranks
570 // (transpose of the preceding communication pattern operation)
571
572 // Build map from local index to post_indices position. Set to -1 for
573 // data that was already on this rank and was therefore was not
574 // sent/received via a postoffice.
575 const std::array<std::int64_t, 2> postoffice_range
576 = dolfinx::MPI::local_range(rank, shape[0], size);
577 std::vector<std::int32_t> post_indices_map(
578 postoffice_range[1] - postoffice_range[0], -1);
579 for (std::size_t i = 0; i < post_indices.size(); ++i)
580 {
581 assert(post_indices[i] < (int)post_indices_map.size());
582 post_indices_map[post_indices[i]] = i;
583 }
584
585 // Build send buffer
586 std::vector<T> send_buffer_data(shape[1] * recv_disp.back());
587 for (std::size_t p = 0; p < recv_disp.size() - 1; ++p)
588 {
589 int offset = recv_disp[p];
590 for (std::int32_t i = recv_disp[p]; i < recv_disp[p + 1]; ++i)
591 {
592 std::int64_t index = recv_buffer_index[i];
593 if (index >= rank_offset and index < (rank_offset + shape0_local))
594 {
595 // I already had this index before any communication
596 std::int32_t local_index = index - rank_offset;
597 std::copy_n(std::next(x.begin(), shape[1] * local_index), shape[1],
598 std::next(send_buffer_data.begin(), shape[1] * offset));
599 }
600 else
601 {
602 // Take from my 'post bag'
603 auto local_index = index - postoffice_range[0];
604 std::int32_t pos = post_indices_map[local_index];
605 assert(pos != -1);
606 std::copy_n(std::next(post_x.begin(), shape[1] * pos), shape[1],
607 std::next(send_buffer_data.begin(), shape[1] * offset));
608 }
609
610 ++offset;
611 }
612 }
613
614 err = MPI_Dist_graph_create_adjacent(
615 comm, src.size(), src.data(), MPI_UNWEIGHTED, dest.size(), dest.data(),
616 MPI_UNWEIGHTED, MPI_INFO_NULL, false, &neigh_comm0);
617 dolfinx::MPI::check_error(comm, err);
618
619 MPI_Datatype compound_type0;
620 MPI_Type_contiguous(shape[1], dolfinx::MPI::mpi_t<T>, &compound_type0);
621 MPI_Type_commit(&compound_type0);
622
623 std::vector<T> recv_buffer_data(shape[1] * send_disp.back());
624 err = MPI_Neighbor_alltoallv(
625 send_buffer_data.data(), num_items_recv.data(), recv_disp.data(),
626 compound_type0, recv_buffer_data.data(), num_items_per_src.data(),
627 send_disp.data(), compound_type0, neigh_comm0);
628 dolfinx::MPI::check_error(comm, err);
629
630 err = MPI_Type_free(&compound_type0);
631 dolfinx::MPI::check_error(comm, err);
632 err = MPI_Comm_free(&neigh_comm0);
633 dolfinx::MPI::check_error(comm, err);
634
635 std::vector<std::int32_t> index_pos_to_buffer(indices.size(), -1);
636 for (std::size_t i = 0; i < src_to_index.size(); ++i)
637 index_pos_to_buffer[std::get<2>(src_to_index[i])] = i;
638
639 // Extra data to return
640 std::vector<T> x_new(shape[1] * indices.size());
641 for (std::size_t i = 0; i < indices.size(); ++i)
642 {
643 const std::int64_t index = indices[i];
644 if (index >= rank_offset and index < (rank_offset + shape0_local))
645 {
646 // Had data from the start in x
647 auto local_index = index - rank_offset;
648 std::copy_n(std::next(x.begin(), shape[1] * local_index), shape[1],
649 std::next(x_new.begin(), shape[1] * i));
650 }
651 else
652 {
653 if (int src = dolfinx::MPI::index_owner(size, index, shape[0]);
654 src == rank)
655 {
656 // In my post office bag
657 auto local_index = index - postoffice_range[0];
658 std::int32_t pos = post_indices_map[local_index];
659 assert(pos != -1);
660 std::copy_n(std::next(post_x.begin(), shape[1] * pos), shape[1],
661 std::next(x_new.begin(), shape[1] * i));
662 }
663 else
664 {
665 // In my received post
666 std::int32_t pos = index_pos_to_buffer[i];
667 assert(pos != -1);
668 std::copy_n(std::next(recv_buffer_data.begin(), shape[1] * pos),
669 shape[1], std::next(x_new.begin(), shape[1] * i));
670 }
671 }
672 }
673
674 return x_new;
675}
676//---------------------------------------------------------------------------
677template <typename U>
678std::vector<typename std::remove_reference_t<typename U::value_type>>
679distribute_data(MPI_Comm comm0, std::span<const std::int64_t> indices,
680 MPI_Comm comm1, const U& x, int shape1)
681{
682 assert(shape1 > 0);
683 assert(x.size() % shape1 == 0);
684 const std::int64_t shape0_local = x.size() / shape1;
685
686 int err;
687 std::int64_t shape0 = 0;
688 err = MPI_Allreduce(&shape0_local, &shape0, 1, MPI_INT64_T, MPI_SUM, comm0);
689 dolfinx::MPI::check_error(comm0, err);
690
691 std::int64_t rank_offset = -1;
692 if (comm1 != MPI_COMM_NULL)
693 {
694 rank_offset = 0;
695 err = MPI_Exscan(&shape0_local, &rank_offset, 1,
696 dolfinx::MPI::mpi_t<std::int64_t>, MPI_SUM, comm1);
697 dolfinx::MPI::check_error(comm1, err);
698 }
699 else
700 {
701 rank_offset = -1;
702 if (!x.empty())
703 throw std::runtime_error("Non-empty data on null MPI communicator");
704 }
705
706 return distribute_from_postoffice(comm0, indices, x, {shape0, shape1},
707 rank_offset);
708}
709//---------------------------------------------------------------------------
710
711} // namespace dolfinx::MPI
Comm(MPI_Comm comm, bool duplicate=true)
Duplicate communicator and wrap duplicate.
Definition MPI.cpp:12
~Comm()
Destructor (frees wrapped communicator)
Definition MPI.cpp:36
MPI_Comm comm() const noexcept
Return the underlying MPI_Comm object.
Definition MPI.cpp:62
Timer for measuring and logging elapsed time durations.
Definition Timer.h:40
MPI support functionality.
Definition MPI.h:30
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:280
std::vector< typename std::remove_reference_t< typename U::value_type > > distribute_data(MPI_Comm comm0, std::span< const std::int64_t > indices, MPI_Comm comm1, const U &x, int shape1)
Distribute rows of a rectangular data array to ranks where they are required (scalable version).
Definition MPI.h:679
std::pair< std::vector< std::int32_t >, std::vector< typename std::remove_reference_t< typename U::value_type > > > distribute_to_postoffice(MPI_Comm comm, const U &x, std::array< std::int64_t, 2 > shape, std::int64_t rank_offset)
Distribute row data to 'post office' ranks.
Definition MPI.h:313
std::vector< int > compute_graph_edges_nbx(MPI_Comm comm, std::span< const int > edges, int tag=static_cast< int >(tag::consensus_nbx))
Determine incoming graph edges using the NBX consensus algorithm.
Definition MPI.cpp:162
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:113
std::vector< typename std::remove_reference_t< typename U::value_type > > distribute_from_postoffice(MPI_Comm comm, std::span< const std::int64_t > indices, const U &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:464
tag
MPI communication tags.
Definition MPI.h:33
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:95
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:80
int size(MPI_Comm comm)
Definition MPI.cpp:72
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
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:89
Definition MPI.h:268
MPI Type.
Definition MPI.h:275