/*! \file ep_gather.cpp \since 2 may 2016 \brief Definitions of MPI collective function: MPI_Gatherv, MPI_Allgatherv */ #include "ep_lib.hpp" #include #include "ep_declaration.hpp" #include "ep_mpi.hpp" using namespace std; namespace ep_lib { 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) { assert(valid_type(datatype)); ::MPI_Aint datasize, lb; ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize); int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first; int num_ep = comm->ep_comm_ptr->size_rank_info[1].second; #pragma omp critical (_gatherv) comm->my_buffer->void_buffer[ep_rank_loc] = const_cast< void* >(sendbuf); MPI_Barrier_local(comm); if(ep_rank_loc == local_root) { for(int i=0; imy_buffer->void_buffer[i], datasize*recvcounts[i]); } MPI_Barrier_local(comm); } int MPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int input_recvcounts[], const int input_displs[], MPI_Datatype recvtype, int root, MPI_Comm comm) { if(!comm->is_ep) return ::MPI_Gatherv(const_cast(sendbuf), sendcount, to_mpi_type(sendtype), recvbuf, const_cast(input_recvcounts), const_cast(input_displs), to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm)); if(comm->is_intercomm) return MPI_Gatherv_intercomm(sendbuf, sendcount, sendtype, recvbuf, input_recvcounts, input_displs, recvtype, root, comm); assert(sendtype == recvtype); int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first; int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first; int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first; int ep_size = comm->ep_comm_ptr->size_rank_info[0].second; int num_ep = comm->ep_comm_ptr->size_rank_info[1].second; int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second; int root_mpi_rank = comm->ep_rank_map->at(root).second; int root_ep_loc = comm->ep_rank_map->at(root).first; ::MPI_Aint datasize, lb; ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize); int *recvcounts; int* displs; recvcounts = new int[ep_size]; displs = new int[ep_size]; bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root; bool is_root = ep_rank == root; void* local_recvbuf; std::vectorlocal_recvcounts(num_ep, 0); std::vectorlocal_displs(num_ep, 0); if(is_root) { copy(input_recvcounts, input_recvcounts+ep_size, recvcounts); copy(input_displs, input_displs+ep_size, displs); } MPI_Bcast(recvcounts, ep_size, MPI_INT, root, comm); MPI_Bcast(displs, ep_size, MPI_INT, root, comm); if(mpi_rank == root_mpi_rank) MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), root_ep_loc, comm); else MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), 0, comm); if(is_master) { int local_recvbuf_size = std::accumulate(local_recvcounts.begin(), local_recvcounts.end(), 0); for(int i=1; i mpi_recvcounts(mpi_size, 0); std::vector mpi_displs(mpi_size, 0); if(is_master) { for(int i=0; iep_rank_map->at(i).second]+=recvcounts[i]; } for(int i=1; impi_comm)); } // reorder data if(is_root) { int offset; for(int i=0; iep_rank_map->at(i).first; j++) if(comm->ep_rank_map->at(i).second == comm->ep_rank_map->at(j).second) { extra += recvcounts[j]; k++; } offset = mpi_displs[comm->ep_rank_map->at(i).second] + extra; memcpy(recvbuf+displs[i]*datasize, tmp_recvbuf+offset*datasize, recvcounts[i]*datasize); } } delete[] recvcounts; delete[] displs; if(is_master) { delete[] local_recvbuf; } if(is_root) delete[] tmp_recvbuf; } int MPI_Gatherv_intercomm(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int input_recvcounts[], const int input_displs[], MPI_Datatype recvtype, int root, MPI_Comm comm) { printf("MPI_Gatherv_intercomm not yet implemented\n"); MPI_Abort(comm, 0); } }