source: XIOS/dev/dev_trunk_omp/extern/ep_dev/ep_scatter.cpp @ 1604

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

branch_openmp merged with trunk r1597

File size: 4.6 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) return ::MPI_Scatter(sendbuf, sendcount, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm));
45    if(comm->is_intercomm) return MPI_Scatter_intercomm(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
46   
47    assert(sendcount == recvcount);
48
49    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
50    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
51    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
52    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
53    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
54    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
55
56    int root_mpi_rank = comm->ep_rank_map->at(root).second;
57    int root_ep_loc = comm->ep_rank_map->at(root).first;
58
59    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
60    bool is_root = ep_rank == root;
61
62    MPI_Datatype datatype = sendtype;
63    int count = sendcount;
64
65    ::MPI_Aint datasize, lb;
66    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
67   
68    void *tmp_sendbuf;
69    if(is_root) tmp_sendbuf = new void*[ep_size * count * datasize];
70
71    // reorder tmp_sendbuf
72    std::vector<int>local_ranks(num_ep);
73    std::vector<int>ranks(ep_size);
74
75    if(mpi_rank == root_mpi_rank) MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), root_ep_loc, comm);
76    else                          MPI_Gather_local(&ep_rank, 1, MPI_INT, local_ranks.data(), 0, comm);
77
78
79    std::vector<int> recvcounts(mpi_size, 0);
80    std::vector<int> displs(mpi_size, 0);
81
82
83    if(is_master)
84    {
85      for(int i=0; i<ep_size; i++)
86      {
87        recvcounts[comm->ep_rank_map->at(i).second]++;
88      }
89
90      for(int i=1; i<mpi_size; i++)
91        displs[i] = displs[i-1] + recvcounts[i-1];
92
93      ::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));
94    }
95
96
97    if(is_root)
98    for(int i=0; i<ep_size; i++)
99    {
100      memcpy(tmp_sendbuf + i*datasize*count, sendbuf + ranks[i]*datasize*count, count*datasize);
101    }
102
103    // MPI_Scatterv from root to masters
104    void* local_recvbuf;
105    if(is_master) local_recvbuf = new void*[datasize * num_ep * count];
106
107
108    if(is_master)
109    {
110      int local_sendcount = num_ep * count;
111      ::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));
112     
113      if(is_root) for(int i=1; i<mpi_size; i++) displs[i] = displs[i-1] + recvcounts[i-1];
114
115      ::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));
116
117    }
118
119    if(mpi_rank == root_mpi_rank) MPI_Scatter_local(local_recvbuf, count, sendtype, recvbuf, recvcount, recvtype, root_ep_loc, comm);
120    else                          MPI_Scatter_local(local_recvbuf, count, sendtype, recvbuf, recvcount, recvtype, 0, comm);
121
122    if(is_root)   delete[] tmp_sendbuf;
123    if(is_master) delete[] local_recvbuf;
124  }
125
126
127  int MPI_Scatter_intercomm(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
128  {
129    printf("MPI_Scatter_intercomm not yet implemented\n");
130    MPI_Abort(comm, 0);
131  }
132
133}
Note: See TracBrowser for help on using the repository browser.