source: XIOS/dev/branch_openmp/extern/ep_dev/ep_gather.cpp @ 1527

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

save dev

File size: 3.9 KB
Line 
1/*!
2   \file ep_gather.cpp
3   \since 2 may 2016
4
5   \brief Definitions of MPI collective function: MPI_Gather, MPI_Allgather
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_Gather_local(const void *sendbuf, int count, MPI_Datatype datatype, void *recvbuf, int local_root, MPI_Comm comm)
19  {
20    assert(valid_type(datatype));
21
22    ::MPI_Aint datasize, lb;
23    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
24
25    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
26    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
27
28    #pragma omp critical (_gather)
29    comm->my_buffer->void_buffer[ep_rank_loc] = const_cast< void* >(sendbuf);
30
31    MPI_Barrier_local(comm);
32
33    if(ep_rank_loc == local_root)
34    {
35      for(int i=0; i<num_ep; i++)
36        memcpy(recvbuf + datasize * i * count, comm->my_buffer->void_buffer[i], datasize * count);
37    }
38
39    MPI_Barrier_local(comm);
40  }
41
42  int MPI_Gather(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_Gather(const_cast<void*>(sendbuf), sendcount, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype),
45                                         root, to_mpi_comm(comm->mpi_comm));
46    if(comm->is_intercomm) return MPI_Gather_intercomm(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
47
48    assert(sendcount == recvcount && sendtype == recvtype);
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    ::MPI_Aint datasize, lb;
61    ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize);
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    void* local_recvbuf;
67
68    if(is_master)
69    {
70      local_recvbuf = new void*[datasize * num_ep * sendcount];
71    }
72
73    void* tmp_recvbuf;
74    if(is_root) tmp_recvbuf = new void*[datasize * recvcount * ep_size];
75
76
77    if(mpi_rank == root_mpi_rank) MPI_Gather_local(sendbuf, sendcount, sendtype, local_recvbuf, root_ep_loc, comm);
78    else                          MPI_Gather_local(sendbuf, sendcount, sendtype, local_recvbuf, 0, comm);
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]+=sendcount;
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_recvbuf, sendcount*num_ep, to_mpi_type(sendtype), tmp_recvbuf, recvcounts.data(), displs.data(), to_mpi_type(recvtype), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
95    }   
96
97
98    // reorder data
99    if(is_root)
100    {
101      int offset;
102      for(int i=0; i<ep_size; i++)
103      {
104        offset = displs[comm->ep_rank_map->at(i).second] + comm->ep_rank_map->at(i).first * sendcount; 
105        memcpy(recvbuf + i*sendcount*datasize, tmp_recvbuf+offset*datasize, sendcount*datasize);
106
107
108      }
109
110    }
111
112
113    if(is_master)
114    {
115      delete[] local_recvbuf;
116    }
117    if(is_root) delete[] tmp_recvbuf;
118   
119  }
120
121
122  int MPI_Gather_intercomm(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
123  {
124    printf("MPI_Gather_intercomm not yet implemented\n");
125    MPI_Abort(comm, 0);
126  }
127}
Note: See TracBrowser for help on using the repository browser.