source: XIOS/dev/dev_trunk_omp/extern/src_ep_dev/ep_gatherv.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: 8.0 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
14using namespace std;
15
16namespace ep_lib
17{
18
19  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)
20  {
21    assert(valid_type(datatype));
22
23    ::MPI_Aint datasize, lb;
24    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
25
26    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
27    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
28
29
30    #pragma omp critical (_gatherv)
31    comm->my_buffer->void_buffer[ep_rank_loc] = const_cast< void* >(sendbuf);
32
33    MPI_Barrier_local(comm);
34
35    if(ep_rank_loc == local_root)
36    {
37      for(int i=0; i<num_ep; i++)
38        memcpy(recvbuf + datasize*displs[i], comm->my_buffer->void_buffer[i], datasize*recvcounts[i]);
39
40    }
41
42    MPI_Barrier_local(comm);
43  }
44
45  int MPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int input_recvcounts[], const int input_displs[],
46                  MPI_Datatype recvtype, int root, MPI_Comm comm)
47  {
48 
49    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),
50                                          to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm));
51    if(comm->is_intercomm) return MPI_Gatherv_intercomm(sendbuf, sendcount, sendtype, recvbuf, input_recvcounts, input_displs, recvtype, root, comm);
52
53    assert(sendtype == recvtype);
54
55   
56    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
57    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
58    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
59    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
60    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
61    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
62
63    int root_mpi_rank = comm->ep_rank_map->at(root).second;
64    int root_ep_loc = comm->ep_rank_map->at(root).first;
65
66    ::MPI_Aint datasize, lb;
67    ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize);
68
69    int *recvcounts;
70    int* displs;
71
72    recvcounts = new int[ep_size];
73    displs = new int[ep_size];
74
75
76    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
77    bool is_root = ep_rank == root;
78
79    void* local_recvbuf;
80    std::vector<int>local_recvcounts(num_ep, 0);
81    std::vector<int>local_displs(num_ep, 0);
82
83
84    if(is_root)
85    { 
86      copy(input_recvcounts, input_recvcounts+ep_size, recvcounts);
87      copy(input_displs, input_displs+ep_size, displs);
88    }
89
90    MPI_Bcast(recvcounts, ep_size, MPI_INT, root, comm);
91    MPI_Bcast(displs, ep_size, MPI_INT, root, comm);
92
93    if(mpi_rank == root_mpi_rank) MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), root_ep_loc, comm);
94    else                          MPI_Gather_local(&sendcount, 1, MPI_INT, local_recvcounts.data(), 0, comm);
95
96   
97    //if(is_root) printf("rank %d =================local_recvcounts = %d %d\n", ep_rank, local_recvcounts[0], local_recvcounts[1]);
98    //if(!is_root && ep_rank_loc==0 && mpi_rank != root_mpi_rank) printf("rank %d =================local_recvcounts = %d %d\n", ep_rank, local_recvcounts[0], local_recvcounts[1]);
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    //if(is_root) printf("rank %d =================local_recvbuf = %d %d\n", ep_rank, static_cast<int*>(local_recvbuf)[0], static_cast<int*>(local_recvbuf)[1]);
115    //if(!is_root && ep_rank_loc==0 && mpi_rank != root_mpi_rank && ep_rank!=6) printf("rank %d =================local_recvbuf = %d %d\n", ep_rank, static_cast<int*>(local_recvbuf)[0], static_cast<int*>(local_recvbuf)[1]);
116    //if(!is_root && ep_rank_loc==0 && mpi_rank != root_mpi_rank && ep_rank==6) printf("rank %d =================local_recvbuf = %d %d %d\n", ep_rank, static_cast<int*>(local_recvbuf)[0], static_cast<int*>(local_recvbuf)[1], static_cast<int*>(local_recvbuf)[2]);
117
118    void* tmp_recvbuf;
119    int tmp_recvbuf_size = std::accumulate(recvcounts, recvcounts+ep_size, 0);
120
121
122
123    if(is_root) tmp_recvbuf = new void*[datasize * tmp_recvbuf_size];
124
125
126    std::vector<int> mpi_recvcounts(mpi_size, 0);
127    std::vector<int> mpi_displs(mpi_size, 0);
128
129
130    if(is_master)
131    {
132      int sendcount_mpi = 0;
133      for(int i=0; i<num_ep; i++)
134      {
135        sendcount_mpi += local_recvcounts[i];
136      }
137
138      for(int i=0; i<ep_size; i++)
139      {
140        mpi_recvcounts[comm->ep_rank_map->at(i).second]+=recvcounts[i];
141      }
142
143      for(int i=1; i<mpi_size; i++)
144        mpi_displs[i] = mpi_displs[i-1] + mpi_recvcounts[i-1];
145
146
147      ::MPI_Gatherv(local_recvbuf, sendcount_mpi, 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));
148      //printf("****************** rank %d, sendcount*num_ep = %d, sendcount_mpi = %d\n", ep_rank, sendcount*num_ep, sendcount_mpi);
149
150      /*if(is_root) printf("rank %d =================tmp_recvbuf = %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d\n", ep_rank, static_cast<int*>(tmp_recvbuf)[0], static_cast<int*>(tmp_recvbuf)[1], static_cast<int*>(tmp_recvbuf)[2],
151                                                                                                                                 static_cast<int*>(tmp_recvbuf)[3], static_cast<int*>(tmp_recvbuf)[4], static_cast<int*>(tmp_recvbuf)[5],
152                                                                                                                                 static_cast<int*>(tmp_recvbuf)[6], static_cast<int*>(tmp_recvbuf)[7], static_cast<int*>(tmp_recvbuf)[8],
153                                                                                                                                 static_cast<int*>(tmp_recvbuf)[9], static_cast<int*>(tmp_recvbuf)[10], static_cast<int*>(tmp_recvbuf)[11],
154                                                                                                                                 static_cast<int*>(tmp_recvbuf)[12], static_cast<int*>(tmp_recvbuf)[13], static_cast<int*>(tmp_recvbuf)[14],
155                                                                                                                                 static_cast<int*>(tmp_recvbuf)[15], static_cast<int*>(tmp_recvbuf)[16]);*/
156    }   
157
158
159    // reorder data
160    if(is_root)
161    {
162      int offset;
163      for(int i=0; i<ep_size; i++)
164      {
165        int extra = 0;
166        for(int j=0, k=0; j<ep_size, k<comm->ep_rank_map->at(i).first; j++)
167          if(comm->ep_rank_map->at(i).second == comm->ep_rank_map->at(j).second)
168          {
169            extra += recvcounts[j];
170            k++;
171          } 
172
173        offset = mpi_displs[comm->ep_rank_map->at(i).second] +  extra;
174
175        memcpy(recvbuf+displs[i]*datasize, tmp_recvbuf+offset*datasize, recvcounts[i]*datasize);
176       
177      }
178
179    }
180
181    delete[] recvcounts;
182    delete[] displs;
183
184    if(is_master)
185    {
186      delete[] local_recvbuf;
187    }
188    if(is_root) delete[] tmp_recvbuf;
189  }
190
191  int MPI_Gatherv_intercomm(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int input_recvcounts[], const int input_displs[],
192                            MPI_Datatype recvtype, int root, MPI_Comm comm)
193  {
194    printf("MPI_Gatherv_intercomm not yet implemented\n");
195    MPI_Abort(comm, 0);
196  }
197
198}
199#endif
Note: See TracBrowser for help on using the repository browser.