source: XIOS/dev/dev_trunk_omp/extern/src_ep_dev/ep_allgatherv.cpp @ 1646

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

branch merged with trunk @1645. arch file (ep&mpi) added for ADA

File size: 3.8 KB
Line 
1#ifdef _usingEP
2/*!
3   \file ep_gather.cpp
4   \since 2 may 2016
5
6   \brief Definitions of MPI collective function: MPI_Gatherv, MPI_Allgatherv
7 */
8
9#include "ep_lib.hpp"
10#include <mpi.h>
11#include "ep_declaration.hpp"
12#include "ep_mpi.hpp"
13
14
15using namespace std;
16
17namespace ep_lib
18{
19
20  int MPI_Allgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[], const int displs[], MPI_Datatype recvtype, MPI_Comm comm)
21  {
22
23    if(!comm->is_ep) return ::MPI_Allgatherv(sendbuf, sendcount, to_mpi_type(sendtype), recvbuf, recvcounts, displs, to_mpi_type(recvtype), to_mpi_comm(comm->mpi_comm));
24    if(comm->is_intercomm) return MPI_Allgatherv_intercomm(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
25
26
27    assert(valid_type(sendtype) && valid_type(recvtype));
28
29    MPI_Datatype datatype = sendtype;
30    int count = sendcount;
31
32    ::MPI_Aint datasize, lb;
33
34    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
35
36
37    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
38    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
39    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
40    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
41    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
42    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
43
44    assert(sendcount == recvcounts[ep_rank]);
45
46    bool is_master = ep_rank_loc==0;
47
48    void* local_recvbuf;
49    void* tmp_recvbuf;
50
51    int recvbuf_size = 0;
52    for(int i=0; i<ep_size; i++)
53      recvbuf_size = max(recvbuf_size, displs[i]+recvcounts[i]);
54
55
56    vector<int>local_recvcounts(num_ep, 0);
57    vector<int>local_displs(num_ep, 0);
58
59    MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), 0, comm);
60    for(int i=1; i<num_ep; i++) local_displs[i] = local_displs[i-1] + local_recvcounts[i-1]; 
61
62
63    if(is_master)
64    {
65      local_recvbuf = new void*[datasize * std::accumulate(local_recvcounts.begin(), local_recvcounts.end(), 0)];
66      tmp_recvbuf = new void*[datasize * std::accumulate(recvcounts, recvcounts+ep_size, 0)];
67    }
68
69    MPI_Gatherv_local(sendbuf, count, datatype, local_recvbuf, local_recvcounts.data(), local_displs.data(), 0, comm);
70
71
72    if(is_master)
73    {
74      std::vector<int>mpi_recvcounts(mpi_size, 0);
75      std::vector<int>mpi_displs(mpi_size, 0);
76
77      int local_sendcount = std::accumulate(local_recvcounts.begin(), local_recvcounts.end(), 0);
78      ::MPI_Allgather(&local_sendcount, 1, to_mpi_type(MPI_INT), mpi_recvcounts.data(), 1, to_mpi_type(MPI_INT), to_mpi_comm(comm->mpi_comm));
79
80      for(int i=1; i<mpi_size; i++)
81        mpi_displs[i] = mpi_displs[i-1] + mpi_recvcounts[i-1];
82
83
84      ::MPI_Allgatherv(local_recvbuf, local_sendcount, to_mpi_type(datatype), tmp_recvbuf, mpi_recvcounts.data(), mpi_displs.data(), to_mpi_type(datatype), to_mpi_comm(comm->mpi_comm));
85
86      // reorder
87      int offset;
88      for(int i=0; i<ep_size; i++)
89      {
90        int extra = 0;
91        for(int j=0, k=0; j<ep_size, k<comm->ep_rank_map->at(i).first; j++)
92          if(comm->ep_rank_map->at(i).second == comm->ep_rank_map->at(j).second)
93          {
94            extra += recvcounts[j];
95            k++;
96          } 
97
98        offset = mpi_displs[comm->ep_rank_map->at(i).second] +  extra;
99
100        memcpy(recvbuf+displs[i]*datasize, tmp_recvbuf+offset*datasize, recvcounts[i]*datasize);
101       
102      }
103
104    }
105
106    MPI_Bcast_local(recvbuf, recvbuf_size, datatype, 0, comm);
107
108    if(is_master)
109    {
110      delete[] local_recvbuf;
111      delete[] tmp_recvbuf;
112    }
113  }
114
115
116  int MPI_Allgatherv_intercomm(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[], const int displs[], MPI_Datatype recvtype, MPI_Comm comm)
117  {
118    printf("MPI_Allgatherv_intercomm not yet implemented\n");
119    MPI_Abort(comm, 0);
120  }
121
122
123}
124#endif
Note: See TracBrowser for help on using the repository browser.