source: XIOS/dev/branch_openmp/extern/src_ep_dev/ep_gatherv.cpp @ 2022

Last change on this file since 2022 was 1642, checked in by yushan, 5 years ago

dev on ADA. add flag switch _usingEP/_usingMPI

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