TRIQS/nda 1.3.0
Multi-dimensional array library for C++
Loading...
Searching...
No Matches
gather.hpp
Go to the documentation of this file.
1// Copyright (c) 2020-2023 Simons Foundation
2//
3// Licensed under the Apache License, Version 2.0 (the "License");
4// you may not use this file except in compliance with the License.
5// You may obtain a copy of the License at
6//
7// http://www.apache.org/licenses/LICENSE-2.0.txt
8//
9// Unless required by applicable law or agreed to in writing, software
10// distributed under the License is distributed on an "AS IS" BASIS,
11// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12// See the License for the specific language governing permissions and
13// limitations under the License.
14//
15// Authors: Olivier Parcollet, Nils Wentzell
16
17/**
18 * @file
19 * @brief Provides an MPI gather function for nda::Array types.
20 */
21
22#pragma once
23
24#include "../basic_functions.hpp"
25#include "../concepts.hpp"
26#include "../exceptions.hpp"
27#include "../traits.hpp"
28
29#include <mpi/mpi.hpp>
30
31#include <type_traits>
32#include <utility>
33#include <vector>
34
35/**
36 * @ingroup av_mpi
37 * @brief Specialization of the `mpi::lazy` class for nda::Array types and the `mpi::tag::gather` tag.
38 *
39 * @details An object of this class is returned when gathering nda::Array objects across multiple MPI processes.
40 *
41 * It models an nda::ArrayInitializer, that means it can be used to initialize and assign to nda::basic_array and
42 * nda::basic_array_view objects. The target array will be a concatenation of the input arrays along the first
43 * dimension (see nda::concatenate).
44 *
45 * See nda::mpi_gather for an example.
46 *
47 * @tparam A nda::Array type to be gathered.
48 */
49template <nda::Array A>
50struct mpi::lazy<mpi::tag::gather, A> {
51 /// Value type of the array/view.
52 using value_type = typename std::decay_t<A>::value_type;
53
54 /// Const view type of the array/view stored in the lazy object.
55 using const_view_type = decltype(std::declval<const A>()());
56
57 /// View of the array/view to be gathered.
58 const_view_type rhs;
59
60 /// MPI communicator.
61 mpi::communicator comm;
62
63 /// MPI root process.
64 const int root{0}; // NOLINT (const is fine here)
65
66 /// Should all processes receive the result.
67 const bool all{false}; // NOLINT (const is fine here)
68
69 /**
70 * @brief Compute the shape of the target array.
71 *
72 * @details It is assumed that the shape of the input array is the same for all MPI processes except for the first
73 * dimension. The target shape will then be the same as the input shape, except that the extent of its first dimension
74 * will be the sum of the extents of the input arrays along the first dimension.
75 *
76 * @warning This makes an MPI call.
77 *
78 * @return Shape of the target array.
79 */
80 [[nodiscard]] auto shape() const {
81 auto dims = rhs.shape();
82 long dim0 = dims[0];
83 if (!all) {
84 dims[0] = mpi::reduce(dim0, comm, root);
85 if (comm.rank() != root) dims[0] = 1;
86 } else
87 dims[0] = mpi::all_reduce(dim0, comm);
88 return dims;
89 }
90
91 /**
92 * @brief Execute the lazy MPI operation and write the result to a target array/view.
93 *
94 * @tparam T nda::Array type of the target array/view.
95 * @param target Target array/view.
96 */
97 template <nda::Array T>
98 void invoke(T &&target) const { // NOLINT (temporary views are allowed here)
99 // check if the arrays can be used in the MPI call
100 if (not target.is_contiguous()) NDA_RUNTIME_ERROR << "Error in MPI gather for nda::Array: Target array needs to be contiguous";
101 static_assert(std::decay_t<A>::layout_t::stride_order_encoded == std::decay_t<T>::layout_t::stride_order_encoded,
102 "Error in MPI gather for nda::Array: Incompatible stride orders");
103
104 // special case for non-mpi runs
105 if (not mpi::has_env) {
106 target = rhs;
107 return;
108 }
109
110 // get target shape and resize or check the target array
111 auto dims = shape();
112 if (all || (comm.rank() == root)) nda::resize_or_check_if_view(target, dims);
113
114 // gather receive counts and memory displacements
115 auto recvcounts = std::vector<int>(comm.size());
116 auto displs = std::vector<int>(comm.size() + 1, 0);
117 int sendcount = rhs.size();
118 auto mpi_int_type = mpi::mpi_type<int>::get();
119 if (!all)
120 MPI_Gather(&sendcount, 1, mpi_int_type, &recvcounts[0], 1, mpi_int_type, root, comm.get());
121 else
122 MPI_Allgather(&sendcount, 1, mpi_int_type, &recvcounts[0], 1, mpi_int_type, comm.get());
123
124 for (int r = 0; r < comm.size(); ++r) displs[r + 1] = recvcounts[r] + displs[r];
125
126 // gather the data
127 auto mpi_value_type = mpi::mpi_type<value_type>::get();
128 if (!all)
129 MPI_Gatherv((void *)rhs.data(), sendcount, mpi_value_type, target.data(), &recvcounts[0], &displs[0], mpi_value_type, root, comm.get());
130 else
131 MPI_Allgatherv((void *)rhs.data(), sendcount, mpi_value_type, target.data(), &recvcounts[0], &displs[0], mpi_value_type, comm.get());
132 }
133};
134
135namespace nda {
136
137 /**
138 * @ingroup av_mpi
139 * @brief Implementation of an MPI gather for nda::basic_array or nda::basic_array_view types.
140 *
141 * @details Since the returned `mpi::lazy` object models an nda::ArrayInitializer, it can be used to initialize/assign
142 * to nda::basic_array and nda::basic_array_view objects:
143 *
144 * @code{.cpp}
145 * // create an array on all processes
146 * nda::array<int, 2> arr(3, 4);
147 *
148 * // ...
149 * // fill array on each process
150 * // ...
151 *
152 * // gather the array to the root process
153 * nda::array<int, 2> res = mpi::gather(arr);
154 * @endcode
155 *
156 * Here, the array `res` will have a shape of `(3 * comm.size(), 4)`.
157 *
158 * @tparam A nda::basic_array or nda::basic_array_view type.
159 * @param a Array or view to be gathered.
160 * @param comm `mpi::communicator` object.
161 * @param root Rank of the root process.
162 * @param all Should all processes receive the result of the gather.
163 * @return An `mpi::lazy` object modelling an nda::ArrayInitializer.
164 */
165 template <typename A>
166 ArrayInitializer<std::remove_reference_t<A>> auto mpi_gather(A &&a, mpi::communicator comm = {}, int root = 0, bool all = false)
167 requires(is_regular_or_view_v<A>)
168 {
169 if (not a.is_contiguous()) NDA_RUNTIME_ERROR << "Error in MPI gather for nda::Array: Array needs to be contiguous";
170 return mpi::lazy<mpi::tag::gather, A>{std::forward<A>(a), comm, root, all};
171 }
172
173} // namespace nda
#define NDA_RUNTIME_ERROR
ArrayInitializer< std::remove_reference_t< A > > auto mpi_gather(A &&a, mpi::communicator comm={}, int root=0, bool all=false)
Implementation of an MPI gather for nda::basic_array or nda::basic_array_view types.
Definition gather.hpp:166
const int root
MPI root process.
Definition gather.hpp:64
auto shape() const
Compute the shape of the target array.
Definition gather.hpp:80
const bool all
Should all processes receive the result.
Definition gather.hpp:67
const_view_type rhs
View of the array/view to be gathered.
Definition gather.hpp:58
void invoke(T &&target) const
Execute the lazy MPI operation and write the result to a target array/view.
Definition gather.hpp:98
mpi::communicator comm
MPI communicator.
Definition gather.hpp:61