source: XIOS/dev/branch_openmp/extern/src_ep_dev/ep_scatterv.cpp @ 1501

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

unify type : MPI_Datatype MPI_Aint

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