source: XIOS/dev/dev_trunk_omp/extern/src_ep_dev/ep_scatterv.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: 5.5 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_Scatterv
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_Scatterv_local(const void *sendbuf, const int sendcounts[], const int displs[], MPI_Datatype sendtype, void *recvbuf, int recvcount,
20                   MPI_Datatype recvtype, int local_root, MPI_Comm comm)
21  {
22
23    assert(valid_type(sendtype) && valid_type(recvtype));
24
25    ::MPI_Aint datasize, lb;
26    ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize);
27
28    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
29    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
30
31    assert(recvcount == sendcounts[ep_rank_loc]);
32
33    if(ep_rank_loc == local_root)
34      comm->my_buffer->void_buffer[local_root] = const_cast<void*>(sendbuf);
35
36    MPI_Barrier_local(comm);
37
38    #pragma omp critical (_scatterv)     
39    memcpy(recvbuf, comm->my_buffer->void_buffer[local_root]+datasize*displs[ep_rank_loc], datasize * recvcount);
40   
41
42    MPI_Barrier_local(comm);
43  }
44
45  int MPI_Scatterv(const void *sendbuf, const int sendcounts[], const int displs[], MPI_Datatype sendtype, void *recvbuf, int recvcount,
46                   MPI_Datatype recvtype, int root, MPI_Comm comm)
47  {
48    if(!comm->is_ep) return ::MPI_Scatterv(sendbuf, sendcounts, displs, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm));
49    if(comm->is_intercomm) return MPI_Scatterv_intercomm(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
50   
51    assert(sendtype == recvtype);
52
53    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
54    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
55    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
56    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
57    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
58    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
59
60    int root_mpi_rank = comm->ep_rank_map->at(root).second;
61    int root_ep_loc = comm->ep_rank_map->at(root).first;
62
63    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
64    bool is_root = ep_rank == root;
65
66    MPI_Datatype datatype = sendtype;
67    int count = recvcount;
68
69    ::MPI_Aint datasize, lb;
70    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
71   
72    void *tmp_sendbuf;
73    if(is_root) tmp_sendbuf = new void*[ep_size * count * datasize];
74
75    // reorder tmp_sendbuf
76    std::vector<int>local_ranks(num_ep);
77    std::vector<int>ranks(ep_size);
78
79    if(mpi_rank == root_mpi_rank) MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), root_ep_loc, comm);
80    else                          MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), 0, comm);
81
82
83    std::vector<int> recvcounts(mpi_size, 0);
84    std::vector<int> my_displs(mpi_size, 0);
85
86
87    if(is_master)
88    {
89      for(int i=0; i<ep_size; i++)
90      {
91        recvcounts[comm->ep_rank_map->at(i).second]++;
92      }
93
94      for(int i=1; i<mpi_size; i++)
95        my_displs[i] = my_displs[i-1] + recvcounts[i-1];
96
97      ::MPI_Gatherv(local_ranks.data(), num_ep, to_mpi_type(MPI_INT), ranks.data(), recvcounts.data(), my_displs.data(), to_mpi_type(MPI_INT), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
98    }
99
100
101    if(is_root)
102    {
103      int local_displs = 0;
104      for(int i=0; i<ep_size; i++)
105      {
106        memcpy(tmp_sendbuf+local_displs, sendbuf + displs[ranks[i]]*datasize, sendcounts[ranks[i]]*datasize);
107        local_displs += sendcounts[ranks[i]]*datasize;
108      }
109    }
110
111    // MPI_Scatterv from root to masters
112 
113    void* local_sendbuf;
114    int local_sendcount;
115    if(mpi_rank == root_mpi_rank) MPI_Reduce_local(&recvcount, &local_sendcount, 1, MPI_INT, MPI_SUM, root_ep_loc, comm);
116    else                          MPI_Reduce_local(&recvcount, &local_sendcount, 1, MPI_INT, MPI_SUM, 0, comm);
117
118    if(is_master) 
119    {
120      local_sendbuf = new void*[datasize * local_sendcount];
121 
122      ::MPI_Gather(&local_sendcount, 1, to_mpi_type(MPI_INT), recvcounts.data(), 1, to_mpi_type(MPI_INT), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
123
124      if(is_root) for(int i=1; i<mpi_size; i++) my_displs[i] = my_displs[i-1] + recvcounts[i-1];
125
126      ::MPI_Scatterv(tmp_sendbuf, recvcounts.data(), my_displs.data(), to_mpi_type(sendtype), local_sendbuf, num_ep*count, to_mpi_type(recvtype), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
127    }
128
129    std::vector<int>local_sendcounts(num_ep, 0);
130    std::vector<int>local_displs(num_ep, 0);
131
132    MPI_Gather_local(&recvcount, 1, MPI_INT, local_sendcounts.data(), 0, comm);
133    MPI_Bcast_local(local_sendcounts.data(), num_ep, MPI_INT, 0, comm);
134    for(int i=1; i<num_ep; i++)
135      local_displs[i] = local_displs[i-1] + local_sendcounts[i-1];
136   
137
138    if(mpi_rank == root_mpi_rank) MPI_Scatterv_local(local_sendbuf, local_sendcounts.data(), local_displs.data(), sendtype, recvbuf, recvcount, recvtype, root_ep_loc, comm);
139    else                          MPI_Scatterv_local(local_sendbuf, local_sendcounts.data(), local_displs.data(), sendtype, recvbuf, recvcount, recvtype, 0, comm);
140
141
142    if(is_root) delete[] tmp_sendbuf;
143    if(is_master) delete[] local_sendbuf;
144  }
145
146
147  int MPI_Scatterv_intercomm(const void *sendbuf, const int sendcounts[], const int displs[], MPI_Datatype sendtype, void *recvbuf, int recvcount,
148                   MPI_Datatype recvtype, int root, MPI_Comm comm)
149  {
150    printf("MPI_Scatterv_intercomm not yet implemented\n");
151    MPI_Abort(comm, 0);
152  }
153
154
155}
156#endif
Note: See TracBrowser for help on using the repository browser.