TRIQS/mpi 1.3.0
C++ interface to MPI
Loading...
Searching...
No Matches
vector.hpp
Go to the documentation of this file.
1// Copyright (c) 2019-2024 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: Thomas Hahn, Olivier Parcollet, Nils Wentzell
16
21
22#pragma once
23
24#include "./communicator.hpp"
26#include "./ranges.hpp"
27#include "./utils.hpp"
28
29#include <mpi.h>
30
31#include <type_traits>
32#include <utility>
33#include <vector>
34
35namespace mpi {
36
41
53 template <typename T> void mpi_broadcast(std::vector<T> &v, communicator c = {}, int root = 0) {
54 auto count = v.size();
55 broadcast(count, c, root);
56 if (c.rank() != root) v.resize(count);
57 broadcast_range(v, c, root);
58 }
59
77 template <typename T> auto mpi_reduce(std::vector<T> const &v, communicator c = {}, int root = 0, bool all = false, MPI_Op op = MPI_SUM) {
78 using value_type = std::remove_cvref_t<decltype(reduce(std::declval<T>()))>;
79 std::vector<value_type> res(c.rank() == root || all ? v.size() : 0);
80 reduce_range(v, res, c, root, all, op);
81 return res;
82 }
83
99 template <typename T1, typename T2>
100 void mpi_reduce_into(std::vector<T1> const &v_in, std::vector<T2> &v_out, communicator c = {}, int root = 0, bool all = false,
101 MPI_Op op = MPI_SUM) {
102 if ((c.rank() == root || all) && v_out.size() != v_in.size()) v_out.resize(v_in.size());
103 reduce_range(v_in, v_out, c, root, all, op);
104 }
105
119 template <typename T> void mpi_scatter_into(std::vector<T> const &v_in, std::vector<T> &v_out, communicator c = {}, int root = 0) {
120 auto scatter_size = static_cast<int>(v_in.size());
121 broadcast(scatter_size, c, root);
122 auto const recvcount = chunk_length(scatter_size, c.size(), c.rank());
123 if (v_out.size() != recvcount) v_out.resize(recvcount);
124 scatter_range(v_in, v_out, scatter_size, c, root);
125 }
126
141 template <typename T> void mpi_gather_into(std::vector<T> const &v_in, std::vector<T> &v_out, communicator c = {}, int root = 0, bool all = false) {
142 auto const gather_size = mpi::all_reduce(v_in.size(), c);
143 if ((c.rank() == root || all) && v_out.size() != gather_size) v_out.resize(gather_size);
144 gather_range(v_in, v_out, c, root, all);
145 }
146
148
149} // namespace mpi
C++ wrapper around MPI_Comm providing various convenience functions.
Provides a C++ wrapper class for an MPI_Comm object.
Provides generic implementations for a subset of collective MPI communications (broadcast,...
auto mpi_reduce(std::array< T, N > const &arr, communicator c={}, int root=0, bool all=false, MPI_Op op=MPI_SUM)
Implementation of an MPI reduce for a std::array.
Definition array.hpp:74
void mpi_reduce_into(std::array< T1, N1 > const &arr_in, std::array< T2, N2 > &arr_out, communicator c={}, int root=0, bool all=false, MPI_Op op=MPI_SUM)
Implementation of an MPI reduce for a std::array that reduces directly into an existing output array.
Definition array.hpp:99
void mpi_scatter_into(std::vector< T > const &v_in, std::vector< T > &v_out, communicator c={}, int root=0)
Implementation of an MPI scatter for a std::vector that scatters directly into an existing output vec...
Definition vector.hpp:119
void mpi_gather_into(T const &x, R &&rg, communicator c={}, int root=0, bool all=false)
Implementation of an MPI gather that gathers directly into an existing output range for types that ha...
void reduce_range(R1 &&in_rg, R2 &&out_rg, communicator c={}, int root=0, bool all=false, MPI_Op op=MPI_SUM)
Implementation of an MPI reduce for std::ranges::sized_range objects.
Definition ranges.hpp:119
void gather_range(R1 &&in_rg, R2 &&out_rg, communicator c={}, int root=0, bool all=false)
Implementation of an MPI gather for mpi::MPICompatibleRange objects.
Definition ranges.hpp:278
void broadcast_range(R &&rg, communicator c={}, int root=0)
Implementation of an MPI broadcast for std::ranges::sized_range objects.
Definition ranges.hpp:69
void scatter_range(R1 &&in_rg, R2 &&out_rg, long scatter_size, communicator c={}, int root=0, long chunk_size=1)
Implementation of an MPI scatter for mpi::MPICompatibleRange objects.
Definition ranges.hpp:213
void mpi_broadcast(std::array< T, N > &arr, communicator c={}, int root=0)
Implementation of an MPI broadcast for a std::array.
Definition array.hpp:53
decltype(auto) all_reduce(T &&x, communicator c={}, MPI_Op op=MPI_SUM)
Generic MPI all-reduce.
void broadcast(T &&x, communicator c={}, int root=0)
Generic MPI broadcast.
long chunk_length(long end, int nranges, int i, long min_size=1)
Get the length of the ith subrange after splitting the integer range [0, end) as evenly as possible a...
Definition chunk.hpp:50
Provides an MPI broadcast, reduce, scatter and gather for generic ranges.
Provides general utilities related to MPI.