/*! \file ep_send.hpp \since 2 may 2016 \brief Definitions of MPI send functions: MPI_Send, MPI_Ssend, MPI_Isend, MPI_Issend */ #include "ep_lib.hpp" #include #include "ep_declaration.hpp" #include "ep_mpi.hpp" namespace ep_lib { int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) { if(!comm->is_ep) return MPI_Send_mpi(buf, count, datatype, dest, tag, comm); if(comm->is_intercomm) dest = comm->inter_rank_map->at(dest); Debug("\nMPI_Send with EP\n"); int ep_src_loc = comm->ep_comm_ptr->size_rank_info[1].first; int ep_dest_loc = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).first; int mpi_tag = tag_combine(tag, ep_src_loc, ep_dest_loc); int mpi_dest = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).second; #ifdef _check_sum check_sum_send(buf, count, datatype, dest, tag, comm); #endif return ::MPI_Send(buf, count, to_mpi_type(datatype), mpi_dest, mpi_tag, to_mpi_comm(comm->mpi_comm)); } int MPI_Ssend(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) { if(!comm->is_ep) return MPI_Ssend_mpi(buf, count, datatype, dest, tag, comm); if(comm->is_intercomm) dest = comm->inter_rank_map->at(dest); Debug("\nMPI_Ssend with EP\n"); int ep_src_loc = comm->ep_comm_ptr->size_rank_info[1].first; int ep_dest_loc = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).first; int mpi_tag = tag_combine(tag, ep_src_loc, ep_dest_loc); int mpi_dest = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).second; #ifdef _check_sum check_sum_send(buf, count, datatype, dest, tag, comm); #endif return ::MPI_Ssend(buf, count, to_mpi_type(datatype), mpi_dest, mpi_tag, to_mpi_comm(comm->mpi_comm)); } int MPI_Isend(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request) { if(!comm->is_ep) return MPI_Isend_mpi(buf, count, datatype, dest, tag, comm, request); if(comm->is_intercomm) dest = comm->inter_rank_map->at(dest); Debug("\nMPI_Isend with EP\n"); int src_rank; MPI_Comm_rank(comm, &src_rank); #ifdef _check_sum check_sum_send(buf, count, datatype, dest, tag, comm); #endif *request = new ep_request; memcheck("new "<< *request <<" : in ep_lib::MPI_Isend, *request = new ep_request"); (*request)->mpi_request = new ::MPI_Request; memcheck("new "<< (*request)->mpi_request <<" : in ep_lib::MPI_Isend, (*request)->mpi_request = new ::MPI_Request"); int ep_src_loc = comm->ep_comm_ptr->size_rank_info[1].first; int ep_dest_loc = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).first; int mpi_tag = tag_combine(tag, ep_src_loc, ep_dest_loc); int mpi_dest = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).second; (*request)->ep_src = src_rank; (*request)->ep_tag = tag; (*request)->ep_datatype = datatype; (*request)->type = 1; // used in wait (*request)->comm = comm; (*request)->buf = const_cast(buf); return ::MPI_Isend(buf, count, to_mpi_type(datatype), mpi_dest, mpi_tag, to_mpi_comm(comm->mpi_comm), to_mpi_request_ptr(*request)); } int MPI_Issend(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request) { if(!comm->is_ep) return MPI_Issend_mpi(buf, count, datatype, dest, tag, comm, request); if(comm->is_intercomm) dest = comm->inter_rank_map->at(dest); Debug("\nMPI_Issend with EP\n"); int src_rank; MPI_Comm_rank(comm, &src_rank); #ifdef _check_sum check_sum_send(buf, count, datatype, dest, tag, comm); #endif *request = new ep_request; memcheck("new "<< *request <<" : in ep_lib::MPI_Issend, *request = new ep_request"); (*request)->mpi_request = new ::MPI_Request; memcheck("new "<< (*request)->mpi_request <<" : in ep_lib::MPI_Issend, (*request)->mpi_request = new ::MPI_Request"); int ep_src_loc = comm->ep_comm_ptr->size_rank_info[1].first; int ep_dest_loc = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).first; int mpi_tag = tag_combine(tag, ep_src_loc, ep_dest_loc); int mpi_dest = comm->ep_comm_ptr->comm_list[0]->ep_rank_map->at(dest).second; (*request)->ep_src = src_rank; (*request)->ep_tag = tag; (*request)->ep_datatype = datatype; (*request)->type = 1; // used in wait (*request)->comm = comm; (*request)->buf = const_cast(buf); return ::MPI_Issend(buf, count, to_mpi_type(datatype), mpi_dest, mpi_tag, to_mpi_comm(comm->mpi_comm), to_mpi_request_ptr(*request)); } int MPI_Send_mpi(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) { Debug("\nMPI_Send with MPI\n"); return ::MPI_Send(buf, count, to_mpi_type(datatype), dest, tag, to_mpi_comm(comm->mpi_comm)); } int MPI_Ssend_mpi(void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) { Debug("\nMPI_Ssend with MPI\n"); return ::MPI_Ssend(buf, count, to_mpi_type(datatype), dest, tag, to_mpi_comm(comm->mpi_comm)); } int MPI_Isend_mpi(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request) { Debug("\nMPI_Isend with MPI\n"); int src_rank; MPI_Comm_rank(comm, &src_rank); *request = new ep_request; memcheck("new "<< *request <<" : in ep_lib::MPI_Isend_mpi, *request = new ep_request"); (*request)->mpi_request = new ::MPI_Request; memcheck("new "<< (*request)->mpi_request <<" : in ep_lib::MPI_Isend_mpi, (*request)->mpi_request = new ::MPI_Request"); (*request)->ep_src = src_rank; (*request)->ep_tag = tag; (*request)->ep_datatype = datatype; (*request)->type = 1; (*request)->comm = comm; return ::MPI_Isend(buf, count, to_mpi_type(datatype), dest, tag, to_mpi_comm(comm->mpi_comm), to_mpi_request_ptr(*request)); } int MPI_Issend_mpi(const void *buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request *request) { Debug("\nMPI_Issend with MPI\n"); int src_rank; MPI_Comm_rank(comm, &src_rank); *request = new ep_request; memcheck("new "<< *request <<" : in ep_lib::MPI_Issend_mpi, *request = new ep_request"); (*request)->mpi_request = new ::MPI_Request; memcheck("new "<< (*request)->mpi_request <<" : in ep_lib::MPI_Issend_mpi, (*request)->mpi_request = new ::MPI_Request"); (*request)->ep_src = src_rank; (*request)->ep_tag = tag; (*request)->ep_datatype = datatype; (*request)->type = 1; (*request)->comm = comm; return ::MPI_Issend(buf, count, to_mpi_type(datatype), dest, tag, to_mpi_comm(comm->mpi_comm), to_mpi_request_ptr(*request)); } }