source: XIOS/dev/branch_openmp/extern/src_ep_dev/ep_gather.cpp @ 1287

Last change on this file since 1287 was 1287, checked in by yushan, 7 years ago

EP updated

File size: 4.0 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      //printf("local_recvbuf = %d %d \n", static_cast<int*>(recvbuf)[0], static_cast<int*>(recvbuf)[1] );
39    }
40
41    MPI_Barrier_local(comm);
42  }
43
44  int MPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
45  {
46    if(!comm.is_ep)
47    {
48      return ::MPI_Gather(const_cast<void*>(sendbuf), sendcount, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype),
49                   root, to_mpi_comm(comm.mpi_comm));
50    }
51
52    assert(sendcount == recvcount && sendtype == recvtype);
53
54    int ep_rank = comm.ep_comm_ptr->size_rank_info[0].first;
55    int ep_rank_loc = comm.ep_comm_ptr->size_rank_info[1].first;
56    int mpi_rank = comm.ep_comm_ptr->size_rank_info[2].first;
57    int ep_size = comm.ep_comm_ptr->size_rank_info[0].second;
58    int num_ep = comm.ep_comm_ptr->size_rank_info[1].second;
59    int mpi_size = comm.ep_comm_ptr->size_rank_info[2].second;
60
61    int root_mpi_rank = comm.rank_map->at(root).second;
62    int root_ep_loc = comm.rank_map->at(root).first;
63
64    ::MPI_Aint datasize, lb;
65    ::MPI_Type_get_extent(to_mpi_type(sendtype), &lb, &datasize);
66
67    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
68    bool is_root = ep_rank == root;
69
70    void* local_recvbuf;
71
72    if(is_master)
73    {
74      local_recvbuf = new void*[datasize * num_ep * sendcount];
75    }
76
77    void* tmp_recvbuf;
78    if(is_root) tmp_recvbuf = new void*[datasize * recvcount * ep_size];
79
80
81    if(mpi_rank == root_mpi_rank) MPI_Gather_local(sendbuf, sendcount, sendtype, local_recvbuf, root_ep_loc, comm);
82    else                          MPI_Gather_local(sendbuf, sendcount, sendtype, local_recvbuf, 0, comm);
83
84    std::vector<int> recvcounts(mpi_size, 0);
85    std::vector<int> displs(mpi_size, 0);
86
87
88    if(is_master)
89    {
90      for(int i=0; i<ep_size; i++)
91      {
92        recvcounts[comm.rank_map->at(i).second]+=sendcount;
93      }
94
95      for(int i=1; i<mpi_size; i++)
96        displs[i] = displs[i-1] + recvcounts[i-1];
97
98      ::MPI_Gatherv(local_recvbuf, sendcount*num_ep, sendtype, tmp_recvbuf, recvcounts.data(), displs.data(), recvtype, root_mpi_rank, to_mpi_comm(comm.mpi_comm));
99    }   
100
101
102    // reorder data
103    if(is_root)
104    {
105      // printf("tmp_recvbuf = %d %d %d %d %d %d %d %d\n", static_cast<int*>(tmp_recvbuf)[0], static_cast<int*>(tmp_recvbuf)[1],
106      //                                                   static_cast<int*>(tmp_recvbuf)[2], static_cast<int*>(tmp_recvbuf)[3],
107      //                                                   static_cast<int*>(tmp_recvbuf)[4], static_cast<int*>(tmp_recvbuf)[5],
108      //                                                   static_cast<int*>(tmp_recvbuf)[6], static_cast<int*>(tmp_recvbuf)[7] );
109
110      int offset;
111      for(int i=0; i<ep_size; i++)
112      {
113        offset = displs[comm.rank_map->at(i).second] + comm.rank_map->at(i).first * sendcount; 
114        memcpy(recvbuf + i*sendcount*datasize, tmp_recvbuf+offset*datasize, sendcount*datasize);
115
116
117      }
118
119    }
120
121
122    if(is_master)
123    {
124      delete[] local_recvbuf;
125    }
126    if(is_root) delete[] tmp_recvbuf;
127   
128  }
129
130}
Note: See TracBrowser for help on using the repository browser.