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