[1381] | 1 | /*! |
---|
| 2 | \file ep_gather.cpp |
---|
| 3 | \since 2 may 2016 |
---|
| 4 | |
---|
| 5 | \brief Definitions of MPI collective function: MPI_Gatherv, MPI_Allgatherv |
---|
| 6 | */ |
---|
| 7 | |
---|
| 8 | #include "ep_lib.hpp" |
---|
| 9 | #include <mpi.h> |
---|
| 10 | #include "ep_declaration.hpp" |
---|
| 11 | #include "ep_mpi.hpp" |
---|
| 12 | |
---|
| 13 | using namespace std; |
---|
| 14 | |
---|
| 15 | namespace ep_lib |
---|
| 16 | { |
---|
| 17 | |
---|
| 18 | int MPI_Gatherv_local(const void *sendbuf, int count, MPI_Datatype datatype, void *recvbuf, const int recvcounts[], const int displs[], int local_root, MPI_Comm comm) |
---|
| 19 | { |
---|
| 20 | assert(valid_type(datatype)); |
---|
| 21 | |
---|
| 22 | ::MPI_Aint datasize, lb; |
---|
| 23 | ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize); |
---|
| 24 | |
---|
[1500] | 25 | int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first; |
---|
| 26 | int num_ep = comm->ep_comm_ptr->size_rank_info[1].second; |
---|
[1381] | 27 | |
---|
| 28 | |
---|
| 29 | #pragma omp critical (_gatherv) |
---|
[1500] | 30 | comm->my_buffer->void_buffer[ep_rank_loc] = const_cast< void* >(sendbuf); |
---|
[1381] | 31 | |
---|
| 32 | MPI_Barrier_local(comm); |
---|
| 33 | |
---|
| 34 | if(ep_rank_loc == local_root) |
---|
| 35 | { |
---|
| 36 | for(int i=0; i<num_ep; i++) |
---|
[1500] | 37 | memcpy(recvbuf + datasize*displs[i], comm->my_buffer->void_buffer[i], datasize*recvcounts[i]); |
---|
[1381] | 38 | |
---|
| 39 | } |
---|
| 40 | |
---|
| 41 | MPI_Barrier_local(comm); |
---|
| 42 | } |
---|
| 43 | |
---|
| 44 | int MPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int input_recvcounts[], const int input_displs[], |
---|
| 45 | MPI_Datatype recvtype, int root, MPI_Comm comm) |
---|
| 46 | { |
---|
| 47 | |
---|
[1500] | 48 | if(!comm->is_ep) |
---|
[1381] | 49 | { |
---|
| 50 | return ::MPI_Gatherv(const_cast<void*>(sendbuf), sendcount, to_mpi_type(sendtype), recvbuf, const_cast<int*>(input_recvcounts), const_cast<int*>(input_displs), |
---|
[1500] | 51 | to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm)); |
---|
[1381] | 52 | } |
---|
| 53 | |
---|
| 54 | |
---|
| 55 | assert(sendtype == recvtype); |
---|
| 56 | |
---|
| 57 | |
---|
[1500] | 58 | int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first; |
---|
| 59 | int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first; |
---|
| 60 | int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first; |
---|
| 61 | int ep_size = comm->ep_comm_ptr->size_rank_info[0].second; |
---|
| 62 | int num_ep = comm->ep_comm_ptr->size_rank_info[1].second; |
---|
| 63 | int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second; |
---|
[1381] | 64 | |
---|
[1500] | 65 | int root_mpi_rank = comm->rank_map->at(root).second; |
---|
| 66 | int root_ep_loc = comm->rank_map->at(root).first; |
---|
[1381] | 67 | |
---|
| 68 | ::MPI_Aint datasize, lb; |
---|
| 69 | ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize); |
---|
| 70 | |
---|
| 71 | int *recvcounts; |
---|
| 72 | int* displs; |
---|
| 73 | |
---|
| 74 | recvcounts = new int[ep_size]; |
---|
| 75 | displs = new int[ep_size]; |
---|
| 76 | |
---|
| 77 | |
---|
| 78 | bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root; |
---|
| 79 | bool is_root = ep_rank == root; |
---|
| 80 | |
---|
| 81 | void* local_recvbuf; |
---|
| 82 | std::vector<int>local_recvcounts(num_ep, 0); |
---|
| 83 | std::vector<int>local_displs(num_ep, 0); |
---|
| 84 | |
---|
| 85 | |
---|
| 86 | if(is_root) |
---|
| 87 | { |
---|
| 88 | copy(input_recvcounts, input_recvcounts+ep_size, recvcounts); |
---|
| 89 | copy(input_displs, input_displs+ep_size, displs); |
---|
| 90 | } |
---|
| 91 | |
---|
| 92 | MPI_Bcast(recvcounts, ep_size, MPI_INT, root, comm); |
---|
| 93 | MPI_Bcast(displs, ep_size, MPI_INT, root, comm); |
---|
| 94 | |
---|
| 95 | if(mpi_rank == root_mpi_rank) MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), root_ep_loc, comm); |
---|
| 96 | else MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), 0, comm); |
---|
| 97 | |
---|
| 98 | |
---|
| 99 | |
---|
| 100 | if(is_master) |
---|
| 101 | { |
---|
| 102 | int local_recvbuf_size = std::accumulate(local_recvcounts.begin(), local_recvcounts.end(), 0); |
---|
| 103 | |
---|
| 104 | for(int i=1; i<num_ep; i++) |
---|
| 105 | local_displs[i] = local_displs[i-1] + local_recvcounts[i-1]; |
---|
| 106 | |
---|
| 107 | local_recvbuf = new void*[datasize * local_recvbuf_size]; |
---|
| 108 | } |
---|
| 109 | |
---|
| 110 | if(mpi_rank == root_mpi_rank) MPI_Gatherv_local(sendbuf, sendcount, sendtype, local_recvbuf, local_recvcounts.data(), local_displs.data(), root_ep_loc, comm); |
---|
| 111 | else MPI_Gatherv_local(sendbuf, sendcount, sendtype, local_recvbuf, local_recvcounts.data(), local_displs.data(), 0, comm); |
---|
| 112 | |
---|
| 113 | |
---|
| 114 | void* tmp_recvbuf; |
---|
| 115 | int tmp_recvbuf_size = std::accumulate(recvcounts, recvcounts+ep_size, 0); |
---|
| 116 | |
---|
| 117 | if(is_root) tmp_recvbuf = new void*[datasize * tmp_recvbuf_size]; |
---|
| 118 | |
---|
| 119 | |
---|
| 120 | std::vector<int> mpi_recvcounts(mpi_size, 0); |
---|
| 121 | std::vector<int> mpi_displs(mpi_size, 0); |
---|
| 122 | |
---|
| 123 | |
---|
| 124 | if(is_master) |
---|
| 125 | { |
---|
| 126 | for(int i=0; i<ep_size; i++) |
---|
| 127 | { |
---|
[1500] | 128 | mpi_recvcounts[comm->rank_map->at(i).second]+=recvcounts[i]; |
---|
[1381] | 129 | } |
---|
| 130 | |
---|
| 131 | for(int i=1; i<mpi_size; i++) |
---|
| 132 | mpi_displs[i] = mpi_displs[i-1] + mpi_recvcounts[i-1]; |
---|
| 133 | |
---|
| 134 | |
---|
[1500] | 135 | ::MPI_Gatherv(local_recvbuf, sendcount*num_ep, to_mpi_type(sendtype), tmp_recvbuf, mpi_recvcounts.data(), mpi_displs.data(), to_mpi_type(recvtype), root_mpi_rank, to_mpi_comm(comm->mpi_comm)); |
---|
[1381] | 136 | } |
---|
| 137 | |
---|
| 138 | |
---|
| 139 | // reorder data |
---|
| 140 | if(is_root) |
---|
| 141 | { |
---|
| 142 | int offset; |
---|
| 143 | for(int i=0; i<ep_size; i++) |
---|
| 144 | { |
---|
| 145 | int extra = 0; |
---|
[1500] | 146 | for(int j=0, k=0; j<ep_size, k<comm->rank_map->at(i).first; j++) |
---|
| 147 | if(comm->rank_map->at(i).second == comm->rank_map->at(j).second) |
---|
[1381] | 148 | { |
---|
| 149 | extra += recvcounts[j]; |
---|
| 150 | k++; |
---|
| 151 | } |
---|
| 152 | |
---|
[1500] | 153 | offset = mpi_displs[comm->rank_map->at(i).second] + extra; |
---|
[1381] | 154 | |
---|
| 155 | memcpy(recvbuf+displs[i]*datasize, tmp_recvbuf+offset*datasize, recvcounts[i]*datasize); |
---|
| 156 | |
---|
| 157 | } |
---|
| 158 | |
---|
| 159 | } |
---|
| 160 | |
---|
| 161 | delete[] recvcounts; |
---|
| 162 | delete[] displs; |
---|
| 163 | |
---|
| 164 | if(is_master) |
---|
| 165 | { |
---|
| 166 | delete[] local_recvbuf; |
---|
| 167 | } |
---|
| 168 | if(is_root) delete[] tmp_recvbuf; |
---|
| 169 | } |
---|
| 170 | |
---|
| 171 | } |
---|