source: XIOS/dev/branch_openmp/extern/ep_dev/ep_scatter.cpp @ 1500

Last change on this file since 1500 was 1500, checked in by yushan, 6 years ago

save dev

File size: 5.8 KB
Line 
1/*!
2   \file ep_gather.cpp
3   \since 2 may 2016
4
5   \brief Definitions of MPI collective function: MPI_Scatter
6 */
7
8#include "ep_lib.hpp"
9#include <mpi.h>
10#include "ep_declaration.hpp"
11#include "ep_mpi.hpp"
12
13using namespace std;
14
15namespace ep_lib
16{
17
18  int MPI_Scatter_local(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int local_root, MPI_Comm comm)
19  {
20    assert(valid_type(sendtype) && valid_type(recvtype));
21    assert(recvcount == sendcount);
22
23    ::MPI_Aint datasize, lb;
24    ::MPI_Type_get_extent(to_mpi_type(sendtype), &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    if(ep_rank_loc == local_root)
31      comm->my_buffer->void_buffer[local_root] = const_cast<void*>(sendbuf);
32
33    MPI_Barrier_local(comm);
34
35    #pragma omp critical (_scatter)     
36    memcpy(recvbuf, comm->my_buffer->void_buffer[local_root]+datasize*ep_rank_loc*sendcount, datasize * recvcount);
37   
38
39    MPI_Barrier_local(comm);
40  }
41
42  int MPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
43  {
44    if(!comm->is_ep)
45    {
46      return ::MPI_Scatter(sendbuf, sendcount, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm));
47    }
48   
49    assert(sendcount == recvcount);
50
51    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
52    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
53    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
54    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
55    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
56    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
57
58    int root_mpi_rank = comm->rank_map->at(root).second;
59    int root_ep_loc = comm->rank_map->at(root).first;
60
61    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
62    bool is_root = ep_rank == root;
63
64    MPI_Datatype datatype = sendtype;
65    int count = sendcount;
66
67    ::MPI_Aint datasize, lb;
68    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
69   
70    void *tmp_sendbuf;
71    if(is_root) tmp_sendbuf = new void*[ep_size * count * datasize];
72
73    // reorder tmp_sendbuf
74    std::vector<int>local_ranks(num_ep);
75    std::vector<int>ranks(ep_size);
76
77    if(mpi_rank == root_mpi_rank) MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), root_ep_loc, comm);
78    else                          MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), 0, comm);
79
80
81    std::vector<int> recvcounts(mpi_size, 0);
82    std::vector<int> displs(mpi_size, 0);
83
84
85    if(is_master)
86    {
87      for(int i=0; i<ep_size; i++)
88      {
89        recvcounts[comm->rank_map->at(i).second]++;
90      }
91
92      for(int i=1; i<mpi_size; i++)
93        displs[i] = displs[i-1] + recvcounts[i-1];
94
95      ::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));
96    }
97
98
99
100    // if(is_root) printf("\nranks = %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d\n", ranks[0], ranks[1], ranks[2], ranks[3], ranks[4], ranks[5], ranks[6], ranks[7],
101    //                                                                                   ranks[8], ranks[9], ranks[10], ranks[11], ranks[12], ranks[13], ranks[14], ranks[15]);
102
103    if(is_root)
104    for(int i=0; i<ep_size; i++)
105    {
106      memcpy(tmp_sendbuf + i*datasize*count, sendbuf + ranks[i]*datasize*count, count*datasize);
107    }
108
109    // if(is_root) printf("\ntmp_sendbuf = %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d %d\n", static_cast<int*>(tmp_sendbuf)[0], static_cast<int*>(tmp_sendbuf)[2], static_cast<int*>(tmp_sendbuf)[4], static_cast<int*>(tmp_sendbuf)[6],
110    //                                                                           static_cast<int*>(tmp_sendbuf)[8], static_cast<int*>(tmp_sendbuf)[10], static_cast<int*>(tmp_sendbuf)[12], static_cast<int*>(tmp_sendbuf)[14],
111    //                                                                           static_cast<int*>(tmp_sendbuf)[16], static_cast<int*>(tmp_sendbuf)[18], static_cast<int*>(tmp_sendbuf)[20], static_cast<int*>(tmp_sendbuf)[22],
112    //                                                                           static_cast<int*>(tmp_sendbuf)[24], static_cast<int*>(tmp_sendbuf)[26], static_cast<int*>(tmp_sendbuf)[28], static_cast<int*>(tmp_sendbuf)[30] );
113
114
115    // MPI_Scatterv from root to masters
116    void* local_recvbuf;
117    if(is_master) local_recvbuf = new void*[datasize * num_ep * count];
118
119
120    if(is_master)
121    {
122      int local_sendcount = num_ep * count;
123      ::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));
124     
125      if(is_root) for(int i=1; i<mpi_size; i++) displs[i] = displs[i-1] + recvcounts[i-1];
126
127      ::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));
128
129      // printf("local_recvbuf = %d %d %d %d\n", static_cast<int*>(local_recvbuf)[0], static_cast<int*>(local_recvbuf)[1], static_cast<int*>(local_recvbuf)[2], static_cast<int*>(local_recvbuf)[3]);
130                                                          // static_cast<int*>(local_recvbuf)[4], static_cast<int*>(local_recvbuf)[5], static_cast<int*>(local_recvbuf)[6], static_cast<int*>(local_recvbuf)[7]);
131    }
132
133    if(mpi_rank == root_mpi_rank) MPI_Scatter_local(local_recvbuf, count, sendtype, recvbuf, recvcount, recvtype, root_ep_loc, comm);
134    else                          MPI_Scatter_local(local_recvbuf, count, sendtype, recvbuf, recvcount, recvtype, 0, comm);
135
136    if(is_root)   delete[] tmp_sendbuf;
137    if(is_master) delete[] local_recvbuf;
138  }
139
140}
Note: See TracBrowser for help on using the repository browser.