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

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

dev on ADA. add flag switch _usingEP/_usingMPI

File size: 5.5 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
[1520]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)
[1520]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)     
[1520]38    memcpy(recvbuf, comm->my_buffer->void_buffer[local_root]+datasize*displs[ep_rank_loc], datasize * recvcount);
[1295]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  {
[1539]47    if(!comm->is_ep) return ::MPI_Scatterv(sendbuf, sendcounts, displs, to_mpi_type(sendtype), recvbuf, recvcount, to_mpi_type(recvtype), root, to_mpi_comm(comm->mpi_comm));
48    if(comm->is_intercomm) return MPI_Scatterv_intercomm(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
[1295]49   
50    assert(sendtype == recvtype);
[1134]51
[1520]52    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
53    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
54    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
55    int ep_size = comm->ep_comm_ptr->size_rank_info[0].second;
56    int num_ep = comm->ep_comm_ptr->size_rank_info[1].second;
57    int mpi_size = comm->ep_comm_ptr->size_rank_info[2].second;
[1134]58
[1520]59    int root_mpi_rank = comm->ep_rank_map->at(root).second;
60    int root_ep_loc = comm->ep_rank_map->at(root).first;
[1134]61
[1295]62    bool is_master = (ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root;
63    bool is_root = ep_rank == root;
[1134]64
[1295]65    MPI_Datatype datatype = sendtype;
66    int count = recvcount;
[1134]67
[1295]68    ::MPI_Aint datasize, lb;
69    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
70   
71    void *tmp_sendbuf;
72    if(is_root) tmp_sendbuf = new void*[ep_size * count * datasize];
[1134]73
[1295]74    // reorder tmp_sendbuf
75    std::vector<int>local_ranks(num_ep);
76    std::vector<int>ranks(ep_size);
[1134]77
[1642]78    if(mpi_rank == root_mpi_rank) MPI_Gather_local(&ep_rank, 1, EP_INT, local_ranks.data(), root_ep_loc, comm);
79    else                          MPI_Gather_local(&ep_rank, 1, EP_INT, local_ranks.data(), 0, comm);
[1134]80
81
[1295]82    std::vector<int> recvcounts(mpi_size, 0);
83    std::vector<int> my_displs(mpi_size, 0);
[1134]84
85
[1295]86    if(is_master)
[1134]87    {
[1295]88      for(int i=0; i<ep_size; i++)
[1134]89      {
[1520]90        recvcounts[comm->ep_rank_map->at(i).second]++;
[1295]91      }
[1134]92
[1295]93      for(int i=1; i<mpi_size; i++)
94        my_displs[i] = my_displs[i-1] + recvcounts[i-1];
[1134]95
[1642]96      ::MPI_Gatherv(local_ranks.data(), num_ep, to_mpi_type(EP_INT), ranks.data(), recvcounts.data(), my_displs.data(), to_mpi_type(EP_INT), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
[1134]97    }
98
99
[1295]100    if(is_root)
[1134]101    {
[1295]102      int local_displs = 0;
103      for(int i=0; i<ep_size; i++)
[1134]104      {
[1295]105        memcpy(tmp_sendbuf+local_displs, sendbuf + displs[ranks[i]]*datasize, sendcounts[ranks[i]]*datasize);
106        local_displs += sendcounts[ranks[i]]*datasize;
[1134]107      }
108    }
109
[1295]110    // MPI_Scatterv from root to masters
111 
112    void* local_sendbuf;
113    int local_sendcount;
[1642]114    if(mpi_rank == root_mpi_rank) MPI_Reduce_local(&recvcount, &local_sendcount, 1, EP_INT, EP_SUM, root_ep_loc, comm);
115    else                          MPI_Reduce_local(&recvcount, &local_sendcount, 1, EP_INT, EP_SUM, 0, comm);
[1134]116
[1295]117    if(is_master) 
[1134]118    {
[1295]119      local_sendbuf = new void*[datasize * local_sendcount];
120 
[1642]121      ::MPI_Gather(&local_sendcount, 1, to_mpi_type(EP_INT), recvcounts.data(), 1, to_mpi_type(EP_INT), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
[1134]122
[1295]123      if(is_root) for(int i=1; i<mpi_size; i++) my_displs[i] = my_displs[i-1] + recvcounts[i-1];
[1134]124
[1520]125      ::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));
[1295]126    }
[1289]127
[1295]128    std::vector<int>local_sendcounts(num_ep, 0);
129    std::vector<int>local_displs(num_ep, 0);
[1289]130
[1642]131    MPI_Gather_local(&recvcount, 1, EP_INT, local_sendcounts.data(), 0, comm);
132    MPI_Bcast_local(local_sendcounts.data(), num_ep, EP_INT, 0, comm);
[1295]133    for(int i=1; i<num_ep; i++)
134      local_displs[i] = local_displs[i-1] + local_sendcounts[i-1];
[1289]135   
[1134]136
[1295]137    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);
138    else                          MPI_Scatterv_local(local_sendbuf, local_sendcounts.data(), local_displs.data(), sendtype, recvbuf, recvcount, recvtype, 0, comm);
[1134]139
[1289]140
[1295]141    if(is_root) delete[] tmp_sendbuf;
142    if(is_master) delete[] local_sendbuf;
143  }
[1289]144
145
[1539]146  int MPI_Scatterv_intercomm(const void *sendbuf, const int sendcounts[], const int displs[], MPI_Datatype sendtype, void *recvbuf, int recvcount,
147                   MPI_Datatype recvtype, int root, MPI_Comm comm)
148  {
149    printf("MPI_Scatterv_intercomm not yet implemented\n");
150    MPI_Abort(comm, 0);
151  }
152
153
[1134]154}
Note: See TracBrowser for help on using the repository browser.