1
0
mirror of https://github.com/paboyle/Grid.git synced 2025-04-09 21:50:45 +01:00

Policy to control async or sync SendRecv

This commit is contained in:
paboyle 2017-02-07 01:24:54 -05:00
parent 61f82216e2
commit 123c673db7

View File

@ -39,9 +39,13 @@ MPI_Comm CartesianCommunicator::communicator_world;
// Should error check all MPI calls. // Should error check all MPI calls.
void CartesianCommunicator::Init(int *argc, char ***argv) { void CartesianCommunicator::Init(int *argc, char ***argv) {
int flag; int flag;
int provided;
MPI_Initialized(&flag); // needed to coexist with other libs apparently MPI_Initialized(&flag); // needed to coexist with other libs apparently
if ( !flag ) { if ( !flag ) {
MPI_Init(argc,argv); // MPI_Init_thread(argc,argv,MPI_THREAD_SERIALIZED,&provided);
// assert (provided == MPI_THREAD_SERIALIZED);
MPI_Init_thread(argc,argv,MPI_THREAD_MULTIPLE,&provided);
assert (provided == MPI_THREAD_MULTIPLE);
} }
MPI_Comm_dup (MPI_COMM_WORLD,&communicator_world); MPI_Comm_dup (MPI_COMM_WORLD,&communicator_world);
ShmInitGeneric(); ShmInitGeneric();
@ -152,24 +156,34 @@ void CartesianCommunicator::SendToRecvFromBegin(std::vector<CommsRequest_t> &lis
int from, int from,
int bytes) int bytes)
{ {
MPI_Request xrq; int myrank = _processor;
MPI_Request rrq;
int rank = _processor;
int ierr; int ierr;
ierr =MPI_Isend(xmit, bytes, MPI_CHAR,dest,_processor,communicator,&xrq); if ( CommunicatorPolicy == CommunicatorPolicyIsend ) {
ierr|=MPI_Irecv(recv, bytes, MPI_CHAR,from,from,communicator,&rrq); MPI_Request xrq;
MPI_Request rrq;
assert(ierr==0);
list.push_back(xrq); ierr =MPI_Irecv(recv, bytes, MPI_CHAR,from,from,communicator,&rrq);
list.push_back(rrq); ierr|=MPI_Isend(xmit, bytes, MPI_CHAR,dest,_processor,communicator,&xrq);
assert(ierr==0);
list.push_back(xrq);
list.push_back(rrq);
} else {
// Give the CPU to MPI immediately; can use threads to overlap optionally
ierr=MPI_Sendrecv(xmit,bytes,MPI_CHAR,dest,myrank,
recv,bytes,MPI_CHAR,from, from,
communicator,MPI_STATUS_IGNORE);
assert(ierr==0);
}
} }
void CartesianCommunicator::SendToRecvFromComplete(std::vector<CommsRequest_t> &list) void CartesianCommunicator::SendToRecvFromComplete(std::vector<CommsRequest_t> &list)
{ {
int nreq=list.size(); if ( CommunicatorPolicy == CommunicatorPolicyIsend ) {
std::vector<MPI_Status> status(nreq); int nreq=list.size();
int ierr = MPI_Waitall(nreq,&list[0],&status[0]); std::vector<MPI_Status> status(nreq);
assert(ierr==0); int ierr = MPI_Waitall(nreq,&list[0],&status[0]);
assert(ierr==0);
}
} }
void CartesianCommunicator::Barrier(void) void CartesianCommunicator::Barrier(void)