On our cluster, ISSM deadlocked during initialization in FemModel::SyncLocalVectorWithClonesVertices. I looked into the code, and there are a few functions that use the same communication pattern:
- all processes send
- then all processes receive what has been sent
this only works as long as MPI has enough space to buffer all the sends, which depends on the MPI implementation and config. If there is not enough buffer space, the sends block until they are received, but no one ever receives anything because all processes are blocking on send.
It's possible that our cluster should be configured with more buffer space just for general performance, but I still feel this is a bug that should be fixed. I made some changes to the functions using MPI_Isend and MPI_Wait that resolve the issue. I attached a patch at the end (I wanted to upload it as a file, but the forum won't let me). A few comments on my changes:
- it's for version 4.20, but should work in new versions.
- I think it would be possible to work with only blocking sends by carefully arranging which process sends first and which receives first. But this would probably be neither faster nor easier
- I didn't handle the Adjoint MPI case as I haven't worked with that before
- there may be an elegant way to put this pattern into a function so it doesn't have to be repeated that often, but I haven't found one so far
Index: src/c/classes/FemModel.cpp
===================================================================
--- src/c/classes/FemModel.cpp (revision 27035)
+++ src/c/classes/FemModel.cpp (working copy)
@@ -1438,7 +1438,13 @@
#ifdef _HAVE_AD_
IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
#else
- IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size());
+ IssmDouble** send_buffers = xNew<IssmDouble*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank = 0;rank<num_procs;rank++){
+ send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ IssmDouble* recv_buffer = xNew<IssmDouble>(this->vertices->Size());
#endif
for(int rank=0;rank<num_procs;rank++){
if(this->vertices->common_send[rank]){
@@ -1447,24 +1453,30 @@
int master_lid = this->vertices->common_send_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(!vertex->clone);
- buffer[i] = local_vector[vertex->lid];
+ send_buffers[rank][i] = local_vector[vertex->lid];
}
- ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->vertices->common_recv[rank]){
int numids = this->vertices->common_recv[rank];
- ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
int master_lid = this->vertices->common_recv_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(vertex->clone);
- local_vector[vertex->lid] = buffer[i];
+ local_vector[vertex->lid] = recv_buffer[i];
}
}
}
- xDelete<IssmDouble>(buffer);
+ xDelete<IssmDouble>(recv_buffer);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<IssmDouble>(send_buffers[rank]);
+ }
+ xDelete<IssmDouble*>(send_buffers);
+ xDelete<ISSM_MPI_Request>(send_requests);
}/*}}}*/
void FemModel::SyncLocalVectorWithClonesVerticesAdd(IssmDouble* local_vector){/*{{{*/
@@ -1477,7 +1489,13 @@
#ifdef _HAVE_AD_
IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size(),"t"); //only one alloc, "t" is required by adolc
#else
- IssmDouble* buffer = xNew<IssmDouble>(this->vertices->Size());
+ IssmDouble** send_buffers = xNew<IssmDouble*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank = 0;rank<num_procs;rank++){
+ send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ IssmDouble* recv_buffer = xNew<IssmDouble>(this->vertices->Size());
#endif
/*1st: add slaves to master values (reverse of what we usually do)*/
@@ -1488,24 +1506,28 @@
int master_lid = this->vertices->common_recv_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(vertex->clone);
- buffer[i] = local_vector[vertex->lid];
+ send_buffers[rank][i] = local_vector[vertex->lid];
}
- ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->vertices->common_send[rank]){
int numids = this->vertices->common_send[rank];
- ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
int master_lid = this->vertices->common_send_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(!vertex->clone);
- local_vector[vertex->lid] += buffer[i];
+ local_vector[vertex->lid] += recv_buffer[i];
}
}
}
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ }
+
/*Now sync masters across partitions*/
for(int rank=0;rank<num_procs;rank++){
if(this->vertices->common_send[rank]){
@@ -1514,24 +1536,30 @@
int master_lid = this->vertices->common_send_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(!vertex->clone);
- buffer[i] = local_vector[vertex->lid];
+ send_buffers[rank][i] = local_vector[vertex->lid];
}
- ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->vertices->common_recv[rank]){
int numids = this->vertices->common_recv[rank];
- ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
int master_lid = this->vertices->common_recv_ids[rank][i];
Vertex* vertex=xDynamicCast<Vertex*>(this->vertices->GetObjectByOffset(master_lid));
_assert_(vertex->clone);
- local_vector[vertex->lid] = buffer[i];
+ local_vector[vertex->lid] = recv_buffer[i];
}
}
}
- xDelete<IssmDouble>(buffer);
+ xDelete<IssmDouble>(recv_buffer);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<IssmDouble>(send_buffers[rank]);
+ }
+ xDelete<IssmDouble*>(send_buffers);
+ xDelete<ISSM_MPI_Request>(send_requests);
}/*}}}*/
void FemModel::GetLocalVectorWithClonesNodes(IssmDouble** plocal_vector,Vector<IssmDouble> *vector){/*{{{*/
@@ -1560,7 +1588,13 @@
#ifdef _HAVE_AD_
IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size(),"t"); //only one alloc, "t" is required by adolc
#else
- IssmDouble* buffer = xNew<IssmDouble>(this->nodes->Size());
+ IssmDouble** send_buffers = xNew<IssmDouble*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank = 0;rank<num_procs;rank++){
+ send_buffers[rank] = xNew<IssmDouble>(this->vertices->Size());
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ IssmDouble* recv_buffer = xNew<IssmDouble>(this->nodes->Size());
#endif
for(int rank=0;rank<num_procs;rank++){
if(this->nodes->common_send[rank]){
@@ -1569,24 +1603,30 @@
int master_lid = this->nodes->common_send_ids[rank][i];
Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
_assert_(!vertex->clone);
- buffer[i] = local_vector[vertex->lid];
+ send_buffers[rank][i] = local_vector[vertex->lid];
}
- ISSM_MPI_Send(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_buffers[rank],numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->nodes->common_recv[rank]){
int numids = this->nodes->common_recv[rank];
- ISSM_MPI_Recv(buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_buffer,numids,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
int master_lid = this->nodes->common_recv_ids[rank][i];
Node* vertex=xDynamicCast<Node*>(this->nodes->GetObjectByOffset(master_lid));
_assert_(vertex->clone);
- local_vector[vertex->lid] = buffer[i];
+ local_vector[vertex->lid] = recv_buffer[i];
}
}
}
- xDelete<IssmDouble>(buffer);
+ xDelete<IssmDouble>(recv_buffer);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<IssmDouble>(send_buffers[rank]);
+ }
+ xDelete<IssmDouble*>(send_buffers);
+ xDelete<ISSM_MPI_Request>(send_requests);
/*Assign output pointer*/
*plocal_vector = local_vector;
Index: src/c/classes/Nodes.cpp
===================================================================
--- src/c/classes/Nodes.cpp (revision 27035)
+++ src/c/classes/Nodes.cpp (working copy)
@@ -180,29 +180,41 @@
/* Finally, remember that cpus may have skipped some objects, because they were clones. For every
* object that is not a clone, tell them to show their dofs, so that later on, they can get picked
* up by their clones: */
- int maxdofspernode = this->MaxNumDofs(GsetEnum);
- int* truedofs = xNewZeroInit<int>(this->Size()*maxdofspernode); //only one alloc
+ int maxdofspernode = this->MaxNumDofs(setenum);
+ int** send_truedofs = xNew<int*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank=0;rank<num_procs;rank++){
+ send_truedofs[rank] = xNewZeroInit<int>(this->Size()*maxdofspernode);
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ int* recv_truedofs = xNewZeroInit<int>(this->Size()*maxdofspernode);
for(int rank=0;rank<num_procs;rank++){
if(this->common_send[rank]){
int numids = this->common_send[rank];
for(int i=0;i<numids;i++){
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
- node->ShowMasterDofs(&truedofs[i*maxdofspernode+0],setenum);
+ node->ShowMasterDofs(&send_truedofs[rank][i*maxdofspernode+0],setenum);
}
- ISSM_MPI_Send(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_truedofs[rank],numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->common_recv[rank]){
int numids = this->common_recv[rank];
- ISSM_MPI_Recv(truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_truedofs,numids*maxdofspernode,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
- node->UpdateCloneDofs(&truedofs[i*maxdofspernode+0],setenum);
+ node->UpdateCloneDofs(&recv_truedofs[i*maxdofspernode+0],setenum);
}
}
}
- xDelete<int>(truedofs);
+ xDelete<int>(recv_truedofs);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<int>(send_truedofs[rank]);
+ }
+ xDelete<int*>(send_truedofs);
+ xDelete<ISSM_MPI_Request>(send_requests);
/*Update indexingupdateflag*/
for(Object* & object : this->objects){
@@ -260,28 +272,40 @@
}
/* Share pids of masters and update pids of clones*/
- int* truepids = xNew<int>(this->Size()); //only one alloc
+ int** send_truepids = xNew<int*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank=0;rank<num_procs;rank++){
+ send_truepids[rank] = xNew<int>(this->Size());
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ int* recv_truepids = xNew<int>(this->Size());
for(int rank=0;rank<num_procs;rank++){
if(this->common_send[rank]){
int numids = this->common_send[rank];
for(int i=0;i<numids;i++){
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
- truepids[i] = node->pid;
+ send_truepids[rank][i] = node->pid;
}
- ISSM_MPI_Send(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_truepids[rank],numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->common_recv[rank]){
int numids = this->common_recv[rank];
- ISSM_MPI_Recv(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
- node->pid = truepids[i];
+ node->pid = recv_truepids[i];
}
}
}
- xDelete<int>(truepids);
+ xDelete<int>(recv_truepids);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<int>(send_truepids[rank]);
+ }
+ xDelete<int*>(send_truepids);
+ xDelete<ISSM_MPI_Request>(send_requests);
/*4. Distribute G dofs once for all*/
//this->DistributeDofs(GsetEnum);
@@ -485,7 +509,13 @@
#ifdef _HAVE_AD_
IssmDouble* buffer = xNew<IssmDouble>(this->Size()*maxdofspernode,"t"); //only one alloc, "t" is required by adolc
#else
- IssmDouble* buffer = xNew<IssmDouble>(this->Size()*maxdofspernode);
+ IssmDouble** send_buffers = xNew<IssmDouble*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank=0;rank<num_procs;rank++){
+ send_buffers[rank] = xNew<IssmDouble>(this->Size()*maxdofspernode);
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ IssmDouble* recv_buffer = xNew<IssmDouble>(this->Size()*maxdofspernode);
#endif
for(int rank=0;rank<num_procs;rank++){
if(this->common_send[rank]){
@@ -494,23 +524,29 @@
int master_lid = this->common_send_ids[rank][i];
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(master_lid));
_assert_(!node->IsClone());
- for(int j=0;j<node->gsize;j++) buffer[i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]];
+ for(int j=0;j<node->gsize;j++) send_buffers[rank][i*maxdofspernode+j]=local_ug[node->gdoflist_local[j]];
}
- ISSM_MPI_Send(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_buffers[rank],numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->common_recv[rank]){
int numids = this->common_recv[rank];
- ISSM_MPI_Recv(buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_buffer,numids*maxdofspernode,ISSM_MPI_DOUBLE,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
int master_lid = this->common_recv_ids[rank][i];
Node* node=xDynamicCast<Node*>(this->GetObjectByOffset(master_lid));
- for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = buffer[i*maxdofspernode+j];
+ for(int j=0;j<node->gsize;j++) local_ug[node->gdoflist_local[j]] = recv_buffer[i*maxdofspernode+j];
}
}
}
- xDelete<IssmDouble>(buffer);
+ xDelete<IssmDouble>(recv_buffer);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<IssmDouble>(send_buffers[rank]);
+ }
+ xDelete<IssmDouble*>(send_buffers);
+ xDelete<ISSM_MPI_Request>(send_requests);
/*Assign output pointer*/
*plocal_ug = local_ug;
Index: src/c/classes/Vertices.cpp
===================================================================
--- src/c/classes/Vertices.cpp (revision 27035)
+++ src/c/classes/Vertices.cpp (working copy)
@@ -234,28 +234,40 @@
}
/* Share pids of masters and update pids of clones*/
- int* truepids = xNew<int>(this->Size()); //only one alloc
+ int** send_truepids = xNew<int*>(num_procs);
+ ISSM_MPI_Request* send_requests = xNew<ISSM_MPI_Request>(num_procs);
+ for (int rank=0;rank<num_procs;rank++){
+ send_truepids[rank] = xNew<int>(this->Size());
+ send_requests[rank] = MPI_REQUEST_NULL;
+ }
+ int* recv_truepids = xNew<int>(this->Size());
for(int rank=0;rank<num_procs;rank++){
if(this->common_send[rank]){
int numids = this->common_send[rank];
for(int i=0;i<numids;i++){
Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_send_ids[rank][i]));
- truepids[i] = vertex->pid;
+ send_truepids[rank][i] = vertex->pid;
}
- ISSM_MPI_Send(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm());
+ ISSM_MPI_Isend(send_truepids[rank],numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&send_requests[rank]);
}
}
for(int rank=0;rank<num_procs;rank++){
if(this->common_recv[rank]){
int numids = this->common_recv[rank];
- ISSM_MPI_Recv(truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
+ ISSM_MPI_Recv(recv_truepids,numids,ISSM_MPI_INT,rank,0,IssmComm::GetComm(),&status);
for(int i=0;i<numids;i++){
Vertex* vertex=xDynamicCast<Vertex*>(this->GetObjectByOffset(this->common_recv_ids[rank][i]));
- vertex->pid = truepids[i];
+ vertex->pid = recv_truepids[i];
}
}
}
- xDelete<int>(truepids);
+ xDelete<int>(recv_truepids);
+ for(int rank=0;rank<num_procs;rank++){
+ ISSM_MPI_Wait(&send_requests[rank],&status);
+ xDelete<int>(send_truepids[rank]);
+ }
+ xDelete<int*>(send_truepids);
+ xDelete<ISSM_MPI_Request>(send_requests);
}/*}}}*/
int Vertices::NumberOfVertices(){/*{{{*/
return this->numberofvertices;
Index: src/c/toolkits/mpi/issmmpi.cpp
===================================================================
--- src/c/toolkits/mpi/issmmpi.cpp (revision 27035)
+++ src/c/toolkits/mpi/issmmpi.cpp (working copy)
@@ -508,6 +508,48 @@
#endif
return rc;
}/*}}}*/
+int ISSM_MPI_Isend(void *buf, int count, ISSM_MPI_Datatype datatype, int dest, int tag, ISSM_MPI_Comm comm, ISSM_MPI_Request* req){ /*{{{*/
+
+ int rc=0;
+#ifdef _HAVE_MPI_
+#if defined(_HAVE_AMPI_) && !defined(_WRAPPERS_)
+ rc=AMPI_Isend(buf,
+ count,
+ datatype,
+ dest,
+ tag,
+ #if !defined(_HAVE_ADJOINTMPI_) && !defined(_HAVE_MEDIPACK_)
+ AMPI_TO_RECV, // as long as there are no other variants
+ #endif
+ comm,
+ req);
+# else
+ rc=MPI_Isend(buf,
+ count,
+ datatype,
+ dest,
+ tag,
+ comm,
+ req);
+# endif
+#else
+// nothing to be done here
+#endif
+ return rc;
+}/*}}}*/
+int ISSM_MPI_Wait(ISSM_MPI_Request *req, ISSM_MPI_Status *status){/*{{{*/
+ int rc=0;
+#ifdef _HAVE_MPI_
+#if defined(_HAVE_AMPI_) && !defined(_WRAPPERS_)
+ rc=AMPI_Wait(req, status);
+# else
+ rc=MPI_Wait(req, status);
+# endif
+#else
+// nothing to be done here
+#endif
+ return rc;
+}
double ISSM_MPI_Wtime(void){/*{{{*/
#ifdef _HAVE_MPI_
Index: src/c/toolkits/mpi/issmmpi.h
===================================================================
--- src/c/toolkits/mpi/issmmpi.h (revision 27035)
+++ src/c/toolkits/mpi/issmmpi.h (working copy)
@@ -71,11 +71,13 @@
typedef AMPI_Datatype ISSM_MPI_Datatype;
typedef AMPI_Op ISSM_MPI_Op;
typedef AMPI_Status ISSM_MPI_Status;
+ typedef AMPI_Request ISSM_MPI_Request;
#else
typedef MPI_Comm ISSM_MPI_Comm;
typedef MPI_Datatype ISSM_MPI_Datatype;
typedef MPI_Op ISSM_MPI_Op;
typedef MPI_Status ISSM_MPI_Status;
+ typedef MPI_Request ISSM_MPI_Request;
#endif
#if defined(_HAVE_MEDIPACK_) && !defined(_WRAPPERS_)
@@ -207,6 +209,8 @@
int ISSM_MPI_Scatter(void *sendbuf, int sendcnt, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcnt, ISSM_MPI_Datatype recvtype, int root, ISSM_MPI_Comm comm);
int ISSM_MPI_Scatterv(void *sendbuf, int *sendcnts, int *displs, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcnt, ISSM_MPI_Datatype recvtype, int root, ISSM_MPI_Comm comm);
int ISSM_MPI_Send(void *buf, int count, ISSM_MPI_Datatype datatype, int dest, int tag, ISSM_MPI_Comm comm);
+int ISSM_MPI_Isend(void* buf, int count, ISSM_MPI_Datatype datatype, int dest, int tag, ISSM_MPI_Comm comm, ISSM_MPI_Request* req);
+int ISSM_MPI_Wait(ISSM_MPI_Request* req, ISSM_MPI_Status* status);
double ISSM_MPI_Wtime(void);
int ISSM_MPI_Comm_split(ISSM_MPI_Comm comm, int color, int key, ISSM_MPI_Comm *newcomm);
int ISSM_MPI_Intercomm_create(ISSM_MPI_Comm comm,int local_leader,ISSM_MPI_Comm peer_comm, int remote_leader, int tag,ISSM_MPI_Comm *newintercomm);