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