Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
MPI_ERR_COMM was sometimes returned instead of MPI_ERR_RANK in SMPI
[simgrid.git] / src / smpi / smpi_pmpi.c
index a83b0bd..5bd1361 100644 (file)
@@ -182,7 +182,6 @@ int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent
 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
 {
   int retval;
-  MPI_Aint dummy;
 
   smpi_bench_end();
   if (datatype == MPI_DATATYPE_NULL) {
@@ -190,7 +189,8 @@ int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
   } else if (extent == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    retval = smpi_datatype_extent(datatype, &dummy, extent);
+    *extent = smpi_datatype_get_extent(datatype);
+    retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
   return retval;
@@ -394,7 +394,7 @@ int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
                            MPI_Group * newgroup)
 {
-  int retval, i, proc1, proc2, size, size2;
+  int retval, i, proc1, proc2, size;
 
   smpi_bench_end();
   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
@@ -402,9 +402,8 @@ int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
   } else if (newgroup == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    size = smpi_group_size(group1);
-    size2 = smpi_group_size(group2);
-    for (i = 0; i < size2; i++) {
+    size = smpi_group_size(group2);
+    for (i = 0; i < size; i++) {
       proc2 = smpi_group_index(group2, i);
       proc1 = smpi_group_rank(group1, proc2);
       if (proc1 == MPI_UNDEFINED) {
@@ -415,12 +414,13 @@ int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
       *newgroup = MPI_GROUP_EMPTY;
     } else {
       *newgroup = smpi_group_new(size);
-      size2 = smpi_group_size(group1);
-      for (i = 0; i < size2; i++) {
-        proc1 = smpi_group_index(group1, i);
-        proc2 = smpi_group_rank(group2, proc1);
-        if (proc2 != MPI_UNDEFINED) {
-          smpi_group_set_mapping(*newgroup, proc1, i);
+      int j=0;
+      for (i = 0; i < smpi_group_size(group2); i++) {
+        proc2 = smpi_group_index(group2, i);
+        proc1 = smpi_group_rank(group1, proc2);
+        if (proc1 != MPI_UNDEFINED) {
+          smpi_group_set_mapping(*newgroup, proc2, j);
+          j++;
         }
       }
     }
@@ -498,7 +498,7 @@ int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
 
 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
 {
-  int retval, i, size, rank, index;
+  int retval, i, j, newsize, oldsize, index;
 
   smpi_bench_end();
   if (group == MPI_GROUP_NULL) {
@@ -511,22 +511,26 @@ int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
     } else if (n == smpi_group_size(group)) {
       *newgroup = MPI_GROUP_EMPTY;
     } else {
-      size = smpi_group_size(group) - n;
-      *newgroup = smpi_group_new(size);
-      rank = 0;
-      while (rank < size) {
-        for (i = 0; i < n; i++) {
-          if (ranks[i] == rank) {
-            break;
-          }
+      oldsize=smpi_group_size(group);
+      newsize = oldsize - n;
+      *newgroup = smpi_group_new(newsize);
+
+      int* to_exclude=xbt_new(int, smpi_group_size(group));
+      for(i=0; i<oldsize; i++)
+        to_exclude[i]=0;
+      for(i=0; i<n; i++)
+        to_exclude[ranks[i]]=1;
+
+      j=0;
+      for(i=0; i<oldsize; i++){
+        if(to_exclude[i]==0){
+          index = smpi_group_index(group, i);
+          smpi_group_set_mapping(*newgroup, index, j);
+          j++;
         }
-        if (i >= n) {
-          index = smpi_group_index(group, rank);
-          smpi_group_set_mapping(*newgroup, index, rank);
-          
-        }
-        rank++;
       }
+
+      xbt_free(to_exclude);
     }
     smpi_group_use(*newgroup);
     retval = MPI_SUCCESS;
@@ -858,6 +862,24 @@ int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
   return retval;
 }
 
+int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
+  int retval;
+
+    smpi_bench_end();
+    if (request == NULL) {
+      retval = MPI_ERR_ARG;
+    } else if (comm == MPI_COMM_NULL) {
+      retval = MPI_ERR_COMM;
+    } else if (dst == MPI_PROC_NULL) {
+      retval = MPI_SUCCESS;
+    } else {
+      *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
+      retval = MPI_SUCCESS;
+    }
+    smpi_bench_begin();
+    return retval;
+}
+
 int PMPI_Start(MPI_Request * request)
 {
   int retval;
@@ -893,9 +915,10 @@ int PMPI_Request_free(MPI_Request * request)
   int retval;
 
   smpi_bench_end();
-  if (request == MPI_REQUEST_NULL) {
+  if (*request == MPI_REQUEST_NULL) {
     retval = MPI_ERR_ARG;
   } else {
+    if((*request)->flags & PERSISTENT)(*request)->refcount--;
     smpi_mpi_request_free(request);
     retval = MPI_SUCCESS;
   }
@@ -963,7 +986,7 @@ int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
     *request = MPI_REQUEST_NULL;
     retval = MPI_SUCCESS;
   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
-    retval = MPI_ERR_COMM;
+    retval = MPI_ERR_RANK;
   } else if (count < 0) {
     retval = MPI_ERR_COUNT;
   } else if (buf==NULL && count > 0) {
@@ -996,7 +1019,50 @@ int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
   return retval;
 }
 
+int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
+  int retval;
 
+  smpi_bench_end();
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (dst == MPI_PROC_NULL) {
+    *request = MPI_REQUEST_NULL;
+    retval = MPI_SUCCESS;
+  } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
+    retval = MPI_ERR_RANK;
+  } else if (count < 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (buf==NULL && count > 0) {
+    retval = MPI_ERR_COUNT;
+  } else if (datatype == MPI_DATATYPE_NULL){
+    retval = MPI_ERR_TYPE;
+  } else if(tag<0 && tag !=  MPI_ANY_TAG){
+    retval = MPI_ERR_TAG;
+  } else {
+
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
+  TRACE_smpi_computing_out(rank);
+  int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
+  TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
+  TRACE_smpi_send(rank, rank, dst_traced);
+#endif
+
+    *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
+    retval = MPI_SUCCESS;
+
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
+  (*request)->send = 1;
+  TRACE_smpi_computing_in(rank);
+#endif
+  }
+
+  smpi_bench_begin();
+  return retval;
+}
 
 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
              MPI_Comm comm, MPI_Status * status)
@@ -1004,7 +1070,6 @@ int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
   int retval;
 
   smpi_bench_end();
-
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else if (src == MPI_PROC_NULL) {
@@ -1012,7 +1077,7 @@ int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
     status->MPI_SOURCE = MPI_PROC_NULL;
     retval = MPI_SUCCESS;
   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
-    retval = MPI_ERR_COMM;
+    retval = MPI_ERR_RANK;
   } else if (count < 0) {
     retval = MPI_ERR_COUNT;
   } else if (buf==NULL && count > 0) {
@@ -1058,7 +1123,7 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
   } else if (dst == MPI_PROC_NULL) {
     retval = MPI_SUCCESS;
   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
-    retval = MPI_ERR_COMM;
+    retval = MPI_ERR_RANK;
   } else if (count < 0) {
     retval = MPI_ERR_COUNT;
   } else if (buf==NULL && count > 0) {
@@ -1090,6 +1155,50 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
   return retval;
 }
 
+
+
+int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
+  int retval;
+
+   smpi_bench_end();
+
+   if (comm == MPI_COMM_NULL) {
+     retval = MPI_ERR_COMM;
+   } else if (dst == MPI_PROC_NULL) {
+     retval = MPI_SUCCESS;
+   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
+     retval = MPI_ERR_RANK;
+   } else if (count < 0) {
+     retval = MPI_ERR_COUNT;
+   } else if (buf==NULL && count > 0) {
+     retval = MPI_ERR_COUNT;
+   } else if (datatype == MPI_DATATYPE_NULL){
+     retval = MPI_ERR_TYPE;
+   } else if(tag<0 && tag !=  MPI_ANY_TAG){
+     retval = MPI_ERR_TAG;
+   } else {
+
+ #ifdef HAVE_TRACING
+   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
+   TRACE_smpi_computing_out(rank);
+   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
+   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
+   TRACE_smpi_send(rank, rank, dst_traced);
+ #endif
+
+     smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
+     retval = MPI_SUCCESS;
+
+ #ifdef HAVE_TRACING
+   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
+   TRACE_smpi_computing_in(rank);
+ #endif
+   }
+
+   smpi_bench_begin();
+   return retval;}
+
+
 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
                  int dst, int sendtag, void *recvbuf, int recvcount,
                  MPI_Datatype recvtype, int src, int recvtag,
@@ -1110,7 +1219,7 @@ int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
       retval = MPI_SUCCESS;
   }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
       (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
-    retval = MPI_ERR_COMM;
+    retval = MPI_ERR_RANK;
   } else if (sendcount < 0 || recvcount<0) {
       retval = MPI_ERR_COUNT;
   } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
@@ -1126,7 +1235,6 @@ int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
   TRACE_smpi_send(rank, rank, dst_traced);
-  TRACE_smpi_send(rank, src_traced, rank);
 #endif
 
 
@@ -1136,7 +1244,6 @@ int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
 
 #ifdef HAVE_TRACING
   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
-  TRACE_smpi_recv(rank, rank, dst_traced);
   TRACE_smpi_recv(rank, src_traced, rank);
   TRACE_smpi_computing_in(rank);
 #endif
@@ -1152,16 +1259,24 @@ int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
                          MPI_Comm comm, MPI_Status * status)
 {
   //TODO: suboptimal implementation
// void *recvbuf;
+  void *recvbuf;
   int retval;
+  if (datatype == MPI_DATATYPE_NULL) {
+      retval = MPI_ERR_TYPE;
+  } else if (count < 0) {
+      retval = MPI_ERR_COUNT;
+  } else {
+    int size = smpi_datatype_get_extent(datatype) * count;
+    recvbuf = xbt_new(char, size);
+    retval =
+        MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
+                     datatype, src, recvtag, comm, status);
+    if(retval==MPI_SUCCESS){
+        smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
+    }
+    xbt_free(recvbuf);
 
-//  size = smpi_datatype_size(datatype) * count;
-//  recvbuf = xbt_new(char, size);
-  retval =
-      MPI_Sendrecv(buf, count, datatype, dst, sendtag, buf, count,
-                   datatype, src, recvtag, comm, status);
-/*memcpy(buf, recvbuf, size * sizeof(char));
-  xbt_free(recvbuf);*/
+  }
   return retval;
 }
 
@@ -1275,9 +1390,8 @@ int PMPI_Wait(MPI_Request * request, MPI_Status * status)
       : -1;
   TRACE_smpi_computing_out(rank);
 
-  MPI_Group group = smpi_comm_group((*request)->comm);
-  int src_traced = smpi_group_index(group, (*request)->src);
-  int dst_traced = smpi_group_index(group, (*request)->dst);
+  int src_traced = (*request)->src;
+  int dst_traced = (*request)->dst;
   int is_wait_for_receive = (*request)->recv;
   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
 #endif
@@ -1307,30 +1421,15 @@ int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
 #ifdef HAVE_TRACING
   //save requests information for tracing
   int i;
-  xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
-  xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
-  xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
+  int *srcs = xbt_new(int, count);
+  int *dsts = xbt_new(int, count);
+  int *recvs = xbt_new(int, count);
   for (i = 0; i < count; i++) {
     MPI_Request req = requests[i];      //already received requests are no longer valid
     if (req) {
-      int *asrc = xbt_new(int, 1);
-      int *adst = xbt_new(int, 1);
-      int *arecv = xbt_new(int, 1);
-      *asrc = req->src;
-      *adst = req->dst;
-      *arecv = req->recv;
-      xbt_dynar_insert_at(srcs, i, asrc);
-      xbt_dynar_insert_at(dsts, i, adst);
-      xbt_dynar_insert_at(recvs, i, arecv);
-      xbt_free(asrc);
-      xbt_free(adst);
-      xbt_free(arecv);
-    } else {
-      int *t = xbt_new(int, 1);
-      xbt_dynar_insert_at(srcs, i, t);
-      xbt_dynar_insert_at(dsts, i, t);
-      xbt_dynar_insert_at(recvs, i, t);
-      xbt_free(t);
+      srcs[i] = req->src;
+      dsts[i] = req->dst;
+      recvs[i] = req->recv;
     }
   }
   int rank_traced = smpi_process_index();
@@ -1347,22 +1446,18 @@ int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
   }
 #ifdef HAVE_TRACING
   if(*index!=MPI_UNDEFINED){
-    int src_traced, dst_traced, is_wait_for_receive;
-    xbt_dynar_get_cpy(srcs, *index, &src_traced);
-    xbt_dynar_get_cpy(dsts, *index, &dst_traced);
-    xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
+    int src_traced = srcs[*index];
+    int dst_traced = dsts[*index];
+    int is_wait_for_receive = recvs[*index];
     if (is_wait_for_receive) {
       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
     }
     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
-    //clean-up of dynars
-    xbt_dynar_free(&srcs);
-    xbt_dynar_free(&dsts);
-    xbt_dynar_free(&recvs);
+    xbt_free(srcs);
+    xbt_free(dsts);
+    xbt_free(recvs);
   }
   TRACE_smpi_computing_in(rank_traced);
-
-
 #endif
   smpi_bench_begin();
   return retval;
@@ -1375,30 +1470,17 @@ int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
 #ifdef HAVE_TRACING
   //save information from requests
   int i;
-  xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
-  xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
-  xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
+  int *srcs = xbt_new(int, count);
+  int *dsts = xbt_new(int, count);
+  int *recvs = xbt_new(int, count);
+  int valid_count = 0;
   for (i = 0; i < count; i++) {
     MPI_Request req = requests[i];
-    if(req){
-      int *asrc = xbt_new(int, 1);
-      int *adst = xbt_new(int, 1);
-      int *arecv = xbt_new(int, 1);
-      *asrc = req->src;
-      *adst = req->dst;
-      *arecv = req->recv;
-      xbt_dynar_insert_at(srcs, i, asrc);
-      xbt_dynar_insert_at(dsts, i, adst);
-      xbt_dynar_insert_at(recvs, i, arecv);
-      xbt_free(asrc);
-      xbt_free(adst);
-      xbt_free(arecv);
-    }else {
-      int *t = xbt_new(int, 1);
-      xbt_dynar_insert_at(srcs, i, t);
-      xbt_dynar_insert_at(dsts, i, t);
-      xbt_dynar_insert_at(recvs, i, t);
-      xbt_free(t);
+    if(req!=MPI_REQUEST_NULL){
+      srcs[valid_count] = req->src;
+      dsts[valid_count] = req->dst;
+      recvs[valid_count] = req->recv;
+      valid_count++;
     }
   }
   int rank_traced = smpi_process_index();
@@ -1408,20 +1490,18 @@ int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
 #endif
   int retval = smpi_mpi_waitall(count, requests, status);
 #ifdef HAVE_TRACING
-  for (i = 0; i < count; i++) {
-    int src_traced, dst_traced, is_wait_for_receive;
-    xbt_dynar_get_cpy(srcs, i, &src_traced);
-    xbt_dynar_get_cpy(dsts, i, &dst_traced);
-    xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
+  for (i = 0; i < valid_count; i++) {
+    int src_traced = srcs[i];
+    int dst_traced = dsts[i];
+    int is_wait_for_receive = recvs[i];
     if (is_wait_for_receive) {
       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
     }
   }
   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
-  //clean-up of dynars
-  xbt_dynar_free(&srcs);
-  xbt_dynar_free(&dsts);
-  xbt_dynar_free(&recvs);
+  xbt_free(srcs);
+  xbt_free(dsts);
+  xbt_free(recvs);
   TRACE_smpi_computing_in(rank_traced);
 #endif
   smpi_bench_begin();
@@ -1475,7 +1555,7 @@ int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm c
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else {
-    smpi_mpi_bcast(buf, count, datatype, root, comm);
+    mpi_coll_bcast_fun(buf, count, datatype, root, comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1499,7 +1579,7 @@ int PMPI_Barrier(MPI_Comm comm)
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else {
-    smpi_mpi_barrier(comm);
+    mpi_coll_barrier_fun(comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1529,7 +1609,7 @@ int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
              || recvtype == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
   } else {
-    smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+    mpi_coll_gather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                     recvtype, root, comm);
     retval = MPI_SUCCESS;
   }
@@ -1592,8 +1672,8 @@ int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
              || recvtype == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
   } else {
-    smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
-                       recvtype, comm);
+    mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+                           recvtype, comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1623,7 +1703,7 @@ int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
   } else if (recvcounts == NULL || displs == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
+    mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
                         displs, recvtype, comm);
     retval = MPI_SUCCESS;
   }
@@ -1655,7 +1735,7 @@ int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
              || recvtype == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
   } else {
-    smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+    mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
                      recvtype, root, comm);
     retval = MPI_SUCCESS;
   }
@@ -1717,7 +1797,7 @@ int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
+    mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1746,7 +1826,7 @@ int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
   } else if (op == MPI_OP_NULL) {
     retval = MPI_ERR_OP;
   } else {
-    smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
+      mpi_coll_allreduce_fun(sendbuf, recvbuf, count, datatype, op, comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1789,12 +1869,10 @@ int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
 {
-  int retval, i, size, count;
-  int *displs;
-  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
-
+  int retval;
   smpi_bench_end();
 #ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
   TRACE_smpi_computing_out(rank);
   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
 #endif
@@ -1807,19 +1885,9 @@ int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
   } else if (recvcounts == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    /* arbitrarily choose root as rank 0 */
-    /* TODO: faster direct implementation ? */
-    size = smpi_comm_size(comm);
-    count = 0;
-    displs = xbt_new(int, size);
-    for (i = 0; i < size; i++) {
-      count += recvcounts[i];
-      displs[i] = 0;
-    }
-    smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
-    smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
-                      recvcounts[rank], datatype, 0, comm);
-    xbt_free(displs);
+
+    mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
+                       datatype,  op, comm);
     retval = MPI_SUCCESS;
   }
 #ifdef HAVE_TRACING
@@ -1834,7 +1902,7 @@ int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
                  MPI_Comm comm)
 {
-  int retval, size, sendsize;
+  int retval;
 
   smpi_bench_end();
 #ifdef HAVE_TRACING
@@ -1848,24 +1916,7 @@ int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
              || recvtype == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
   } else {
-    size = smpi_comm_size(comm);
-    sendsize = smpi_datatype_size(sendtype) * sendcount;
-    if (sendsize < 200 && size > 12) {
-      retval =
-          smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
-                                         recvbuf, recvcount, recvtype,
-                                         comm);
-    } else if (sendsize < 3000) {
-      retval =
-          smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
-                                                sendtype, recvbuf,
-                                                recvcount, recvtype, comm);
-    } else {
-      retval =
-          smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
-                                            recvbuf, recvcount, recvtype,
-                                            comm);
-    }
+    retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
   }
 #ifdef HAVE_TRACING
   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
@@ -1897,7 +1948,7 @@ int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
     retval = MPI_ERR_ARG;
   } else {
     retval =
-        smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
+        mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
                                   recvbuf, recvcounts, recvdisps, recvtype,
                                   comm);
   }
@@ -1960,7 +2011,7 @@ int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_typ
   } else if (count<0){
     retval = MPI_ERR_COUNT;
   } else {
-    retval = smpi_datatype_contiguous(count, old_type, new_type);
+    retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
   }
   smpi_bench_begin();
   return retval;
@@ -2055,6 +2106,17 @@ int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype*
   smpi_bench_begin();
   return retval;}
 
+int PMPI_Error_class(int errorcode, int* errorclass) {
+  // assume smpi uses only standard mpi error codes
+  *errorclass=errorcode;
+  return MPI_SUCCESS;
+}
+
+
+int PMPI_Initialized(int* flag) {
+   *flag=(smpi_process_data()!=NULL);
+   return MPI_SUCCESS;
+}
 
 /* The following calls are not yet implemented and will fail at runtime. */
 /* Once implemented, please move them above this notice. */
@@ -2128,10 +2190,6 @@ int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
    return not_yet_implemented();
 }
 
-int PMPI_Error_class(int errorcode, int* errorclass) {
-   return not_yet_implemented();
-}
-
 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
    return not_yet_implemented();
 }
@@ -2183,13 +2241,7 @@ int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcou
    return not_yet_implemented();
 }
 
-int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
-   return not_yet_implemented();
-}
 
-int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
-   return not_yet_implemented();
-}
 
 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
    return not_yet_implemented();
@@ -2219,11 +2271,6 @@ int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
    return not_yet_implemented();
 }
 
-int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
-   return not_yet_implemented();
-}
-
-
 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
    return not_yet_implemented();
 }
@@ -2272,6 +2319,32 @@ int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
    return not_yet_implemented();
 }
 
-int PMPI_Initialized(int* flag) {
+int PMPI_Win_fence( int assert,  MPI_Win win){
    return not_yet_implemented();
 }
+
+int PMPI_Win_free( MPI_Win* win){
+   return not_yet_implemented();
+}
+
+int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
+  return not_yet_implemented();
+}
+
+int PMPI_Info_create( MPI_Info *info){
+  return not_yet_implemented();
+}
+
+int PMPI_Info_set( MPI_Info *info, char *key, char *value){
+  return not_yet_implemented();
+}
+
+int PMPI_Info_free( MPI_Info *info){
+  return not_yet_implemented();
+}
+
+int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
+    MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
+  return not_yet_implemented();
+}
+