Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Kill SMPE: use visualization instead.
[simgrid.git] / src / smpi / smpi_mpi.c
index 2b9caeb..674a1ed 100644 (file)
+/* Copyright (c) 2007, 2008, 2009, 2010. The SimGrid Team.
+ * All rights reserved.                                                     */
 
+/* This program is free software; you can redistribute it and/or modify it
+  * under the terms of the license (GNU LGPL) which comes with this package. */
 
 #include "private.h"
 #include "smpi_coll_private.h"
 #include "smpi_mpi_dt_private.h"
 
+#define MPI_CALL_IMPLEM(type,name,args) \
+   type name args __attribute__((weak,alias("P"#name))); \
+   type P##name args
+
 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
                                 "Logging specific to SMPI (mpi)");
 
-int SMPI_MPI_Init(int *argc, char ***argv)
+#ifdef HAVE_TRACING
+//this function need to be here because of the calls to smpi_bench
+int TRACE_smpi_set_category(const char *category)
+{
+  //need to end bench otherwise categories for execution tasks are wrong
+  smpi_bench_end();
+  int ret;
+  if (!IS_TRACING){
+    ret = 1;
+  }else{
+    if (category != NULL) {
+      ret = TRACE_category(category);
+      TRACE_category_set(SIMIX_process_self(), category);
+    }else{
+      //if category is NULL, trace of platform is disabled for this process
+      TRACE_category_unset(SIMIX_process_self());
+      ret = 0;
+    }
+  }
+  //begin bench after changing process's category
+  smpi_bench_begin();
+  return ret;
+}
+#endif
+
+/* MPI User level calls */
+
+MPI_CALL_IMPLEM(int, MPI_Init, (int *argc, char ***argv))
 {
   smpi_process_init(argc, argv);
+#ifdef HAVE_TRACING
+  TRACE_smpi_init(smpi_process_index());
+#endif
+  smpi_bench_begin();
+  return MPI_SUCCESS;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Finalize, (void))
+{
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  TRACE_smpi_finalize(smpi_process_index());
+#endif
+  smpi_process_destroy();
+  return MPI_SUCCESS;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Init_thread, (int *argc, char ***argv, int required, int *provided))
+{
+  if (provided != NULL) {
+    *provided = MPI_THREAD_MULTIPLE;
+  }
+  return MPI_Init(argc, argv);
+}
+
+MPI_CALL_IMPLEM(int, MPI_Query_thread, (int *provided))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (provided == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *provided = MPI_THREAD_MULTIPLE;
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Is_thread_main, (int *flag))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (flag == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *flag = smpi_process_index() == 0;
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Abort, (MPI_Comm comm, int errorcode))
+{
+  smpi_bench_end();
+  smpi_process_destroy();
+  // FIXME: should kill all processes in comm instead
+  SIMIX_process_kill(SIMIX_process_self());
+  return MPI_SUCCESS;
+}
+
+double MPI_Wtime(void)
+{
+  double time;
+
+  smpi_bench_end();
+  time = SIMIX_get_clock();
+  smpi_bench_begin();
+  return time;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Address, (void *location, MPI_Aint * address))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (!address) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *address = (MPI_Aint) location;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_free, (MPI_Datatype * datatype))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (!datatype) {
+    retval = MPI_ERR_ARG;
+  } else {
+    // FIXME: always fail for now
+    retval = MPI_ERR_TYPE;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_size, (MPI_Datatype datatype, int *size))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (size == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *size = (int) smpi_datatype_size(datatype);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_get_extent,
+                       (MPI_Datatype datatype, MPI_Aint * lb,
+                        MPI_Aint * extent))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (lb == NULL || extent == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    retval = smpi_datatype_extent(datatype, lb, extent);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_extent, (MPI_Datatype datatype, MPI_Aint * extent))
+{
+  int retval;
+  MPI_Aint dummy;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (extent == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    retval = smpi_datatype_extent(datatype, &dummy, extent);
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_lb, (MPI_Datatype datatype, MPI_Aint * disp))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (disp == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *disp = smpi_datatype_lb(datatype);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Type_ub, (MPI_Datatype datatype, MPI_Aint * disp))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (disp == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *disp = smpi_datatype_ub(datatype);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Op_create, (MPI_User_function * function, int commute, MPI_Op * op))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (function == NULL || op == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *op = smpi_op_new(function, commute);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Op_free, (MPI_Op * op))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (op == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (*op == MPI_OP_NULL) {
+    retval = MPI_ERR_OP;
+  } else {
+    smpi_op_destroy(*op);
+    *op = MPI_OP_NULL;
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_free, (MPI_Group * group))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (group == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_group_destroy(*group);
+    *group = MPI_GROUP_NULL;
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_size, (MPI_Group group, int *size))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (size == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *size = smpi_group_size(group);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_rank, (MPI_Group group, int *rank))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (rank == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *rank = smpi_group_rank(group, smpi_process_index());
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_translate_ranks,
+                       (MPI_Group group1, int n, int *ranks1,
+                        MPI_Group group2, int *ranks2))
+{
+  int retval, i, index;
+
+  smpi_bench_end();
+  if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else {
+    for (i = 0; i < n; i++) {
+      index = smpi_group_index(group1, ranks1[i]);
+      ranks2[i] = smpi_group_rank(group2, index);
+    }
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_compare, (MPI_Group group1, MPI_Group group2, int *result))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (result == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *result = smpi_group_compare(group1, group2);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_union,
+                       (MPI_Group group1, MPI_Group group2,
+                        MPI_Group * newgroup))
+{
+  int retval, i, proc1, proc2, size, size2;
+
+  smpi_bench_end();
+  if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } 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++) {
+      proc2 = smpi_group_index(group2, i);
+      proc1 = smpi_group_rank(group1, proc2);
+      if (proc1 == MPI_UNDEFINED) {
+        size++;
+      }
+    }
+    if (size == 0) {
+      *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);
+        smpi_group_set_mapping(*newgroup, proc1, i);
+      }
+      for (i = size2; i < size; i++) {
+        proc2 = smpi_group_index(group2, i - size2);
+        smpi_group_set_mapping(*newgroup, proc2, i);
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_intersection,
+                       (MPI_Group group1, MPI_Group group2,
+                        MPI_Group * newgroup))
+{
+  int retval, i, proc1, proc2, size, size2;
+
+  smpi_bench_end();
+  if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } 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++) {
+      proc2 = smpi_group_index(group2, i);
+      proc1 = smpi_group_rank(group1, proc2);
+      if (proc1 == MPI_UNDEFINED) {
+        size--;
+      }
+    }
+    if (size == 0) {
+      *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);
+        }
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_difference,
+                       (MPI_Group group1, MPI_Group group2,
+                        MPI_Group * newgroup))
+{
+  int retval, i, proc1, proc2, size, size2;
+
+  smpi_bench_end();
+  if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newgroup == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    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) {
+        size--;
+      }
+    }
+    if (size == 0) {
+      *newgroup = MPI_GROUP_EMPTY;
+    } else {
+      *newgroup = smpi_group_new(size);
+      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);
+        }
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_incl,
+                       (MPI_Group group, int n, int *ranks,
+                        MPI_Group * newgroup))
+{
+  int retval, i, index;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newgroup == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    if (n == 0) {
+      *newgroup = MPI_GROUP_EMPTY;
+    } else if (n == smpi_group_size(group)) {
+      *newgroup = group;
+    } else {
+      *newgroup = smpi_group_new(n);
+      for (i = 0; i < n; i++) {
+        index = smpi_group_index(group, ranks[i]);
+        smpi_group_set_mapping(*newgroup, index, i);
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_excl,
+                       (MPI_Group group, int n, int *ranks,
+                        MPI_Group * newgroup))
+{
+  int retval, i, size, rank, index;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newgroup == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    if (n == 0) {
+      *newgroup = group;
+    } 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;
+          }
+        }
+        if (i >= n) {
+          index = smpi_group_index(group, rank);
+          smpi_group_set_mapping(*newgroup, index, rank);
+          rank++;
+        }
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_range_incl,
+                       (MPI_Group group, int n, int ranges[][3],
+                        MPI_Group * newgroup))
+{
+  int retval, i, j, rank, size, index;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newgroup == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    if (n == 0) {
+      *newgroup = MPI_GROUP_EMPTY;
+    } else {
+      size = 0;
+      for (i = 0; i < n; i++) {
+        for (rank = ranges[i][0];       /* First */
+             rank >= 0 && rank <= ranges[i][1]; /* Last */
+             rank += ranges[i][2] /* Stride */ ) {
+          size++;
+        }
+      }
+      if (size == smpi_group_size(group)) {
+        *newgroup = group;
+      } else {
+        *newgroup = smpi_group_new(size);
+        j = 0;
+        for (i = 0; i < n; i++) {
+          for (rank = ranges[i][0];     /* First */
+               rank >= 0 && rank <= ranges[i][1];       /* Last */
+               rank += ranges[i][2] /* Stride */ ) {
+            index = smpi_group_index(group, rank);
+            smpi_group_set_mapping(*newgroup, index, j);
+            j++;
+          }
+        }
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Group_range_excl,
+                       (MPI_Group group, int n, int ranges[][3],
+                        MPI_Group * newgroup))
+{
+  int retval, i, newrank, rank, size, index, add;
+
+  smpi_bench_end();
+  if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newgroup == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    if (n == 0) {
+      *newgroup = group;
+    } else {
+      size = smpi_group_size(group);
+      for (i = 0; i < n; i++) {
+        for (rank = ranges[i][0];       /* First */
+             rank >= 0 && rank <= ranges[i][1]; /* Last */
+             rank += ranges[i][2] /* Stride */ ) {
+          size--;
+        }
+      }
+      if (size == 0) {
+        *newgroup = MPI_GROUP_EMPTY;
+      } else {
+        *newgroup = smpi_group_new(size);
+        newrank = 0;
+        while (newrank < size) {
+          for (i = 0; i < n; i++) {
+            add = 1;
+            for (rank = ranges[i][0];   /* First */
+                 rank >= 0 && rank <= ranges[i][1];     /* Last */
+                 rank += ranges[i][2] /* Stride */ ) {
+              if (rank == newrank) {
+                add = 0;
+                break;
+              }
+            }
+            if (add == 1) {
+              index = smpi_group_index(group, newrank);
+              smpi_group_set_mapping(*newgroup, index, newrank);
+            }
+          }
+        }
+      }
+    }
+    smpi_group_use(*newgroup);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_rank, (MPI_Comm comm, int *rank))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    *rank = smpi_comm_rank(comm);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_size, (MPI_Comm comm, int *size))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (size == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *size = smpi_comm_size(comm);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_group, (MPI_Comm comm, MPI_Group * group))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (group == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *group = smpi_comm_group(comm);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_compare, (MPI_Comm comm1, MPI_Comm comm2, int *result))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (result == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    if (comm1 == comm2) {       /* Same communicators means same groups */
+      *result = MPI_IDENT;
+    } else {
+      *result =
+          smpi_group_compare(smpi_comm_group(comm1),
+                             smpi_comm_group(comm2));
+      if (*result == MPI_IDENT) {
+        *result = MPI_CONGRUENT;
+      }
+    }
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_dup, (MPI_Comm comm, MPI_Comm * newcomm))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (newcomm == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *newcomm = smpi_comm_new(smpi_comm_group(comm));
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_create, (MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (group == MPI_GROUP_NULL) {
+    retval = MPI_ERR_GROUP;
+  } else if (newcomm == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *newcomm = smpi_comm_new(group);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_free, (MPI_Comm * comm))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (*comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    smpi_comm_destroy(*comm);
+    *comm = MPI_COMM_NULL;
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Comm_split, (MPI_Comm comm, int color, int key, MPI_Comm* comm_out))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (comm_out == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    *comm_out = smpi_comm_split(comm, color, key);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Send_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 {
+    *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Recv_init,
+                       (void *buf, int count, MPI_Datatype datatype, int src,
+                        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 {
+    *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Start, (MPI_Request * request))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_start(*request);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Startall, (int count, MPI_Request * requests))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (requests == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_startall(count, requests);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Request_free, (MPI_Request * request))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_request_free(request);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Irecv,
+                       (void *buf, int count, MPI_Datatype datatype, int src,
+                        int tag, MPI_Comm comm, MPI_Request * request))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
+#endif
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
+  (*request)->recv = 1;
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Isend,
+                       (void *buf, int count, MPI_Datatype datatype, int dst,
+                        int tag, MPI_Comm comm, MPI_Request * request))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
+  TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
+  TRACE_smpi_send(rank, rank, dst_traced);
+#endif
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
+  (*request)->send = 1;
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Recv,
+                       (void *buf, int count, MPI_Datatype datatype, int src, int tag,
+                        MPI_Comm comm, MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
+  TRACE_smpi_recv(rank, src_traced, rank);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Send,
+                       (void *buf, int count, MPI_Datatype datatype, int dst, int tag,
+                        MPI_Comm comm))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
+  TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
+  TRACE_smpi_send(rank, rank, dst_traced);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    smpi_mpi_send(buf, count, datatype, dst, tag, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Sendrecv,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        int dst, int sendtag, void *recvbuf, int recvcount,
+                        MPI_Datatype recvtype, int src, int recvtag,
+                        MPI_Comm comm, MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
+  int src_traced = smpi_group_rank(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
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else {
+    smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
+                      recvcount, recvtype, src, recvtag, comm, status);
+    retval = MPI_SUCCESS;
+  }
+#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);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Sendrecv_replace,
+                       (void *buf, int count, MPI_Datatype datatype,
+                        int dst, int sendtag, int src, int recvtag,
+                        MPI_Comm comm, MPI_Status * status))
+{
+  //TODO: suboptimal implementation
+  void *recvbuf;
+  int retval, size;
+
+  size = smpi_datatype_size(datatype) * count;
+  recvbuf = xbt_new(char, size);
+  retval =
+      MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
+                   datatype, src, recvtag, comm, status);
+  memcpy(buf, recvbuf, size * sizeof(char));
+  xbt_free(recvbuf);
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Test, (MPI_Request * request, int *flag, MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (request == NULL || flag == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (*request == MPI_REQUEST_NULL) {
+    retval = MPI_ERR_REQUEST;
+  } else {
+    *flag = smpi_mpi_test(request, status);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Testany,
+                       (int count, MPI_Request requests[], int *index, int *flag,
+                        MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+  if (index == NULL || flag == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *flag = smpi_mpi_testany(count, requests, index, status);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Wait, (MPI_Request * request, MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = request && (*request)->comm != MPI_COMM_NULL
+      ? smpi_comm_rank((*request)->comm)
+      : -1;
+  MPI_Group group = smpi_comm_group((*request)->comm);
+  int src_traced = smpi_group_rank(group, (*request)->src);
+  int dst_traced = smpi_group_rank(group, (*request)->dst);
+  int is_wait_for_receive = (*request)->recv;
+  TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
+#endif
+  if (request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (*request == MPI_REQUEST_NULL) {
+    retval = MPI_ERR_REQUEST;
+  } else {
+    smpi_mpi_wait(request, status);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
+  if (is_wait_for_receive) {
+    TRACE_smpi_recv(rank, src_traced, dst_traced);
+  }
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Waitany,
+                       (int count, MPI_Request requests[], int *index,
+                        MPI_Status * status))
+{
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  //save requests information for tracing
+  int i;
+  xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
+  xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
+  xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
+  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);
+    } 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);
+    }
+  }
+  int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
+  TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
+#endif
+  if (index == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *index = smpi_mpi_waitany(count, requests, status);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  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);
+  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_free(srcs);
+  xbt_free(dsts);
+  xbt_free(recvs);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
+
+MPI_CALL_IMPLEM(int, MPI_Waitall, (int count, MPI_Request requests[], MPI_Status status[]))
+{
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  //save information from requests
+  int i;
+  xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
+  xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
+  xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
+  for (i = 0; i < count; i++) {
+    MPI_Request req = requests[i];      //all req should be valid in Waitall
+    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);
+  }
+  int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
+  TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
+#endif
+  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);
+    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_free(srcs);
+  xbt_free(dsts);
+  xbt_free(recvs);
+#endif
   smpi_bench_begin();
   return MPI_SUCCESS;
 }
 
-int SMPI_MPI_Finalize()
+MPI_CALL_IMPLEM(int, MPI_Waitsome,
+                       (int incount, MPI_Request requests[], int *outcount,
+                        int *indices, MPI_Status status[]))
 {
+  int retval;
+
   smpi_bench_end();
-  smpi_process_finalize();
-  return MPI_SUCCESS;
+  if (outcount == NULL || indices == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
 }
 
-// right now this just exits the current node, should send abort signal to all
-// hosts in the communicator (TODO)
-int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
+MPI_CALL_IMPLEM(int, MPI_Bcast,
+                       (void *buf, int count, MPI_Datatype datatype, int root,
+                        MPI_Comm comm))
 {
-  smpi_exit(errorcode);
-  return 0;
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else {
+    smpi_mpi_bcast(buf, count, datatype, root, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
 }
 
-int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
+MPI_CALL_IMPLEM(int, MPI_Barrier, (MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
+  int retval;
 
   smpi_bench_end();
-
-  if (NULL == comm) {
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
-  } else if (NULL == size) {
-    retval = MPI_ERR_ARG;
   } else {
-    *size = comm->size;
+    smpi_mpi_barrier(comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
+MPI_CALL_IMPLEM(int, MPI_Gather,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
+                        int root, MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
+  int retval;
 
   smpi_bench_end();
-
-  if (NULL == comm) {
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
-  } else if (NULL == rank) {
-    retval = MPI_ERR_ARG;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
   } else {
-    *rank = smpi_mpi_comm_rank(comm);
+    smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+                    recvtype, root, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-
-
-/**
- * Barrier
- **/
-int SMPI_MPI_Barrier(MPI_Comm comm)
+MPI_CALL_IMPLEM(int, MPI_Gatherv,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int *recvcounts, int *displs,
+                        MPI_Datatype recvtype, int root, MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
-  int arity=4;
+  int retval;
 
   smpi_bench_end();
-
-  if (NULL == comm) {
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (recvcounts == NULL || displs == NULL) {
+    retval = MPI_ERR_ARG;
   } else {
-
-    /*
-     * original implemantation:
-     * retval = smpi_mpi_barrier(comm);
-     * this one is unrealistic: it just cond_waits, means no time.
-     */
-     retval = nary_tree_barrier( comm, arity );
+    smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
+                     displs, recvtype, root, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-
-
-int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
-                   int tag, MPI_Comm comm, MPI_Request * request)
+MPI_CALL_IMPLEM(int, MPI_Allgather,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
+                        MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
+  int retval;
 
   smpi_bench_end();
-
-  retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
-                               request);
-  if (NULL != *request && MPI_SUCCESS == retval) {
-    retval = smpi_mpi_irecv(*request);
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else {
+    smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+                       recvtype, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
-                  int tag, MPI_Comm comm, MPI_Status * status)
+MPI_CALL_IMPLEM(int, MPI_Allgatherv,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int *recvcounts, int *displs,
+                        MPI_Datatype recvtype, MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
-  smpi_mpi_request_t request;
+  int retval;
 
   smpi_bench_end();
-
-  retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
-                               &request);
-  if (NULL != request && MPI_SUCCESS == retval) {
-    retval = smpi_mpi_irecv(request);
-    if (MPI_SUCCESS == retval) {
-      retval = smpi_mpi_wait(request, status);
-    }
-    xbt_mallocator_release(smpi_global->request_mallocator, request);
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (recvcounts == NULL || displs == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
+                        displs, recvtype, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
-                   int tag, MPI_Comm comm, MPI_Request * request)
+MPI_CALL_IMPLEM(int, MPI_Scatter,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
+                        int root, MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
+  int retval;
 
   smpi_bench_end();
-
-  retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
-                               request);
-  if (NULL != *request && MPI_SUCCESS == retval) {
-    retval = smpi_mpi_isend(*request);
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else {
+    smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
+                     recvtype, root, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-/**
- * MPI_Send user level
- **/
-int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
-                  int tag, MPI_Comm comm)
+MPI_CALL_IMPLEM(int, MPI_Scatterv,
+                       (void *sendbuf, int *sendcounts, int *displs,
+                        MPI_Datatype sendtype, void *recvbuf, int recvcount,
+                        MPI_Datatype recvtype, int root, MPI_Comm comm))
 {
-  int retval = MPI_SUCCESS;
-  smpi_mpi_request_t request;
+  int retval;
 
   smpi_bench_end();
-
-  retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
-                               &request);
-  if (NULL != request && MPI_SUCCESS == retval) {
-    retval = smpi_mpi_isend(request);
-    if (MPI_SUCCESS == retval) {
-      smpi_mpi_wait(request, MPI_STATUS_IGNORE);
-    }
-    xbt_mallocator_release(smpi_global->request_mallocator, request);
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (sendcounts == NULL || displs == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
+                      recvcount, recvtype, root, comm);
+    retval = MPI_SUCCESS;
   }
-
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
 }
 
-
-/**
- * MPI_Sendrecv internal level 
- **/
-int smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
-                   void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
-                   MPI_Comm comm, MPI_Status *status)
+MPI_CALL_IMPLEM(int, MPI_Reduce,
+                       (void *sendbuf, void *recvbuf, int count,
+                        MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm))
 {
-int rank;
-int retval = MPI_SUCCESS;
-smpi_mpi_request_t srequest;
-smpi_mpi_request_t rrequest;
-
-         rank = smpi_mpi_comm_rank(comm);
-
-         /* send */
-         retval = smpi_create_request(sendbuf, sendcount, sendtype, 
-                               rank,dest,sendtag, 
-                               comm, &srequest);
-         smpi_mpi_isend(srequest);
-        
-         /* recv */
-         retval = smpi_create_request(recvbuf, recvcount, recvtype, 
-                               source, rank,recvtag, 
-                               comm, &rrequest);
-         smpi_mpi_irecv(rrequest);
-
-         smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
-         //printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
-         smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
-         //printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
+  int retval;
 
-         return(retval);
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } 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);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
 }
-/**
- * MPI_Sendrecv user entry point
- **/
-int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
-                   void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
-                   MPI_Comm comm, MPI_Status *status)
+
+MPI_CALL_IMPLEM(int, MPI_Allreduce,
+                       (void *sendbuf, void *recvbuf, int count,
+                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm))
 {
-int retval = MPI_SUCCESS;
+  int retval;
 
   smpi_bench_end();
-  smpi_mpi_sendrecv( sendbuf, sendcount, sendtype, dest, sendtag, 
-                        recvbuf, recvcount, recvtype, source, recvtag,
-                        comm, status);
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (op == MPI_OP_NULL) {
+    retval = MPI_ERR_OP;
+  } else {
+    smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
   smpi_bench_begin();
-
   return retval;
-
-       
 }
 
-/**
- * MPI_Wait and friends
- **/
-int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
+MPI_CALL_IMPLEM(int, MPI_Scan,
+                       (void *sendbuf, void *recvbuf, int count,
+                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm))
 {
-  return smpi_mpi_wait(*request, status);
+  int retval;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (op == MPI_OP_NULL) {
+    retval = MPI_ERR_OP;
+  } else {
+    smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
 }
 
-int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
+MPI_CALL_IMPLEM(int, MPI_Reduce_scatter,
+                       (void *sendbuf, void *recvbuf, int *recvcounts,
+                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm))
 {
-  return smpi_mpi_waitall(count, requests, status);
+  int retval, i, size, count;
+  int *displs;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (op == MPI_OP_NULL) {
+    retval = MPI_ERR_OP;
+  } 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);
+    retval = MPI_SUCCESS;
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
 }
 
-int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
-                     MPI_Status status[])
+MPI_CALL_IMPLEM(int, MPI_Alltoall,
+                       (void *sendbuf, int sendcount, MPI_Datatype sendtype,
+                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
+                        MPI_Comm comm))
 {
-  return smpi_mpi_waitany(count, requests, index, status);
-}
+  int retval, size, sendsize;
 
-/**
- * MPI_Bcast
- **/
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || 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);
+    }
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
 
-/**
- * flat bcast 
- **/
-int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
-int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
-                MPI_Comm comm)
+MPI_CALL_IMPLEM(int, MPI_Alltoallv,
+                       (void *sendbuf, int *sendcounts, int *senddisps,
+                        MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
+                        int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm))
 {
-        int rank;
-        int retval = MPI_SUCCESS;
-        smpi_mpi_request_t request;
+  int retval;
 
-        rank = smpi_mpi_comm_rank(comm);
-        if (rank == root) {
-                retval = smpi_create_request(buf, count, datatype, root,
-                                (root + 1) % comm->size, 0, comm, &request);
-                request->forward = comm->size - 1;
-                smpi_mpi_isend(request);
-        } else {
-                retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
-                                0, comm, &request);
-                smpi_mpi_irecv(request);
-        }
-
-        smpi_mpi_wait(request, MPI_STATUS_IGNORE);
-        xbt_mallocator_release(smpi_global->request_mallocator, request);
+  smpi_bench_end();
+#ifdef HAVE_TRACING
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+  TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
+#endif
+  if (comm == MPI_COMM_NULL) {
+    retval = MPI_ERR_COMM;
+  } else if (sendtype == MPI_DATATYPE_NULL
+             || recvtype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
+             || recvdisps == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    retval =
+        smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
+                                  recvbuf, recvcounts, recvdisps, recvtype,
+                                  comm);
+  }
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+#endif
+  smpi_bench_begin();
+  return retval;
+}
 
-        return(retval);
 
-}
-/**
- * Bcast internal level 
- **/
-int smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
-                   MPI_Comm comm)
+MPI_CALL_IMPLEM(int, MPI_Get_processor_name, (char *name, int *resultlen))
 {
   int retval = MPI_SUCCESS;
-  //retval = flat_tree_bcast(buf, count, datatype, root, comm);
-  retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
+
+  smpi_bench_end();
+  strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
+          MPI_MAX_PROCESSOR_NAME - 1);
+  *resultlen =
+      strlen(name) >
+      MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
+
+  smpi_bench_begin();
   return retval;
 }
 
-/**
- * Bcast user entry point
- **/
-int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
-                   MPI_Comm comm)
+MPI_CALL_IMPLEM(int, MPI_Get_count, (MPI_Status * status, MPI_Datatype datatype, int *count))
 {
   int retval = MPI_SUCCESS;
+  size_t size;
 
   smpi_bench_end();
-  smpi_mpi_bcast(buf,count,datatype,root,comm);
+  if (status == NULL || count == NULL) {
+    retval = MPI_ERR_ARG;
+  } else if (datatype == MPI_DATATYPE_NULL) {
+    retval = MPI_ERR_TYPE;
+  } else {
+    size = smpi_datatype_size(datatype);
+    if (size == 0) {
+      *count = 0;
+    } else if (status->count % size != 0) {
+      retval = MPI_UNDEFINED;
+    } else {
+      *count = smpi_mpi_get_count(status, datatype);
+    }
+  }
   smpi_bench_begin();
-
   return retval;
 }
 
+/* The following calls are not yet implemented and will fail at runtime. */
+/* Once implemented, please move them above this notice. */
 
+static int not_yet_implemented(void) {
+   xbt_die("Not yet implemented");
+   return MPI_ERR_OTHER;
+}
 
-#ifdef DEBUG_REDUCE
-/**
- * debugging helper function
- **/
-static void print_buffer_int(void *buf, int len, char *msg, int rank)
-{
-  int tmp, *v;
-  printf("**[%d] %s: ", rank, msg);
-  for (tmp = 0; tmp < len; tmp++) {
-    v = buf;
-    printf("[%d]", v[tmp]);
-  }
-  printf("\n");
-  free(msg);
+MPI_CALL_IMPLEM(int, MPI_Pack_size, (int incount, MPI_Datatype datatype, MPI_Comm comm, int* size)) {
+   return not_yet_implemented();
 }
-static void print_buffer_double(void *buf, int len, char *msg, int rank)
-{
-  int tmp;
-  double *v;
-  printf("**[%d] %s: ", rank, msg);
-  for (tmp = 0; tmp < len; tmp++) {
-    v = buf;
-    printf("[%lf]", v[tmp]);
-  }
-  printf("\n");
-  free(msg);
+
+MPI_CALL_IMPLEM(int, MPI_Cart_coords, (MPI_Comm comm, int rank, int maxdims, int* coords)) {
+   return not_yet_implemented();
 }
 
+MPI_CALL_IMPLEM(int, MPI_Cart_create, (MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart)) {
+   return not_yet_implemented();
+}
 
-#endif
-/**
- * MPI_Reduce internal level 
- **/
-int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
-                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
-{
-        int retval = MPI_SUCCESS;
-        int rank;
-        int size;
-        int i;
-        int tag = 0;
-        smpi_mpi_request_t *requests;
-        smpi_mpi_request_t request;
+MPI_CALL_IMPLEM(int, MPI_Cart_get, (MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords)) {
+   return not_yet_implemented();
+}
 
-        smpi_bench_end();
+MPI_CALL_IMPLEM(int, MPI_Cart_map, (MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank)) {
+   return not_yet_implemented();
+}
 
-        rank = smpi_mpi_comm_rank(comm);
-        size = comm->size;
+MPI_CALL_IMPLEM(int, MPI_Cart_rank, (MPI_Comm comm, int* coords, int* rank)) {
+   return not_yet_implemented();
+}
 
-        if (rank != root) {           // if i am not ROOT, simply send my buffer to root
+MPI_CALL_IMPLEM(int, MPI_Cart_shift, (MPI_Comm comm, int direction, int displ, int* source, int* dest)) {
+   return not_yet_implemented();
+}
 
-#ifdef DEBUG_REDUCE
-                print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
-#endif
-                retval = smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
-                                        &request);
-                smpi_mpi_isend(request);
-                smpi_mpi_wait(request, MPI_STATUS_IGNORE);
-                xbt_mallocator_release(smpi_global->request_mallocator, request);
+MPI_CALL_IMPLEM(int, MPI_Cart_sub, (MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new)) {
+   return not_yet_implemented();
+}
 
-        } else {
-                // i am the ROOT: wait for all buffers by creating one request by sender
-                int src;
-                requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
+MPI_CALL_IMPLEM(int, MPI_Cartdim_get, (MPI_Comm comm, int* ndims)) {
+   return not_yet_implemented();
+}
 
-                void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
-                for (i = 0; i < size-1; i++) {
-                        // we need 1 buffer per request to store intermediate receptions
-                        tmpbufs[i] = xbt_malloc(count * datatype->size);
-                }  
-                // root: initiliaze recv buf with my own snd buf
-                memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
+MPI_CALL_IMPLEM(int, MPI_Graph_create, (MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph)) {
+   return not_yet_implemented();
+}
 
-                // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
-                // since we should op values as soon as one receiving request matches.
-                for (i = 0; i < size-1; i++) {
-                        // reminder: for smpi_create_request() the src is always the process sending.
-                        src = i < root ? i : i + 1;
-                        retval = smpi_create_request(tmpbufs[i], count, datatype,
-                                        src, root, tag, comm, &(requests[i]));
-                        if (NULL != requests[i] && MPI_SUCCESS == retval) {
-                                if (MPI_SUCCESS == retval) {
-                                        smpi_mpi_irecv(requests[i]);
-                                }
-                        }
-                }
-                // now, wait for completion of all irecv's.
-                for (i = 0; i < size-1; i++) {
-                        int index = MPI_UNDEFINED;
-                        smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
-#ifdef DEBUG_REDUCE
-                        printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
-                        print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
-                                        rank);
-#endif
-
-                        // arg 2 is modified
-                        op->func(tmpbufs[index], recvbuf, &count, &datatype);
-#ifdef DEBUG_REDUCE
-                        print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
-#endif
-                        xbt_free(tmpbufs[index]);
-                        /* FIXME: with the following line, it  generates an
-                         * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
-                         */
-                        // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
-                }
-                xbt_free(requests);
-                xbt_free(tmpbufs);
-        }
-        return retval;
+MPI_CALL_IMPLEM(int, MPI_Graph_get, (MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges)) {
+   return not_yet_implemented();
 }
 
-/**
- * MPI_Reduce user entry point
- **/
-int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
-                   MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
-{
-int retval = MPI_SUCCESS;
+MPI_CALL_IMPLEM(int, MPI_Graph_map, (MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank)) {
+   return not_yet_implemented();
+}
 
-         smpi_bench_end();
+MPI_CALL_IMPLEM(int, MPI_Graph_neighbors, (MPI_Comm comm, int rank, int maxneighbors, int* neighbors)) {
+   return not_yet_implemented();
+}
 
-         retval = smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
+MPI_CALL_IMPLEM(int, MPI_Graph_neighbors_count, (MPI_Comm comm, int rank, int* nneighbors)) {
+   return not_yet_implemented();
+}
 
-         smpi_bench_begin();
-         return retval;
+MPI_CALL_IMPLEM(int, MPI_Graphdims_get, (MPI_Comm comm, int* nnodes, int* nedges)) {
+   return not_yet_implemented();
 }
 
+MPI_CALL_IMPLEM(int, MPI_Topo_test, (MPI_Comm comm, int* top_type)) {
+   return not_yet_implemented();
+}
 
+MPI_CALL_IMPLEM(int, MPI_Error_class, (int errorcode, int* errorclass)) {
+   return not_yet_implemented();
+}
 
-/**
- * MPI_Allreduce
- *
- * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
- **/
-int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
-                   MPI_Op op, MPI_Comm comm )
-{
-int retval = MPI_SUCCESS;
-int root=0;  // arbitrary choice
+MPI_CALL_IMPLEM(int, MPI_Errhandler_create, (MPI_Handler_function* function, MPI_Errhandler* errhandler)) {
+   return not_yet_implemented();
+}
 
-         smpi_bench_end();
+MPI_CALL_IMPLEM(int, MPI_Errhandler_free, (MPI_Errhandler* errhandler)) {
+   return not_yet_implemented();
+}
 
-         retval = smpi_mpi_reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
-         if (MPI_SUCCESS != retval)
-                   return(retval);
+MPI_CALL_IMPLEM(int, MPI_Errhandler_get, (MPI_Comm comm, MPI_Errhandler* errhandler)) {
+   return not_yet_implemented();
+}
 
-         retval = smpi_mpi_bcast( sendbuf, count, datatype, root, comm);
+MPI_CALL_IMPLEM(int, MPI_Error_string, (int errorcode, char* string, int* resultlen)) {
+   return not_yet_implemented();
+}
 
-         smpi_bench_end();
-         return( retval );
+MPI_CALL_IMPLEM(int, MPI_Errhandler_set, (MPI_Comm comm, MPI_Errhandler errhandler)) {
+   return not_yet_implemented();
 }
 
+MPI_CALL_IMPLEM(int, MPI_Type_contiguous, (int count, MPI_Datatype old_type, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
+}
 
-/**
- * MPI_Scatter user entry point
- **/
-int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
-                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
-                          int root, MPI_Comm comm)
-{
-  int retval = MPI_SUCCESS;
-  int i;
-  int cnt=0;  
-  int rank;
-  int tag=0;
-  char *cptr;  // to manipulate the void * buffers
-  smpi_mpi_request_t *requests;
-  smpi_mpi_request_t request;
-  smpi_mpi_status_t status;
-
-
-  smpi_bench_end();
-
-  rank = smpi_mpi_comm_rank(comm);
-
-  requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
-  if (rank == root) {
-          // i am the root: distribute my sendbuf
-          //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
-          cptr = sendbuf;
-          for (i=0; i < comm->size; i++) {
-                  if ( i!=root ) { // send to processes ...
-
-                          retval = smpi_create_request((void *)cptr, sendcount, 
-                                          datatype, root, i, tag, comm, &(requests[cnt]));
-                          if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
-                                  if (MPI_SUCCESS == retval) {
-                                          smpi_mpi_isend(requests[cnt]);
-                                  }
-                                 }
-                                 cnt++;
-                       } 
-                       else { // ... except if it's me.
-                                 memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
-                       }
-                  cptr += sendcount*datatype->size;
-           }
-           for(i=0; i<cnt; i++) { // wait for send to complete
-                           /* FIXME: waitall() should be slightly better */
-                           smpi_mpi_wait(requests[i], &status);
-                           xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
-
-           }
-  } 
-  else {  // i am a non-root process: wait data from the root
-           retval = smpi_create_request(recvbuf,recvcount, 
-                                 recvtype, root, rank, tag, comm, &request);
-           if (NULL != request && MPI_SUCCESS == retval) {
-                       if (MPI_SUCCESS == retval) {
-                                 smpi_mpi_irecv(request);
-                       }
-           }
-           smpi_mpi_wait(request, &status);
-           xbt_mallocator_release(smpi_global->request_mallocator, request);
-  }
-  xbt_free(requests);
-
-  smpi_bench_begin();
-
-  return retval;
-}
-
-
-/**
- * MPI_Alltoall user entry point
- * 
- * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
- * ompi/mca/coll/tuned/coll_tuned_module.c
- **/
-int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
-                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
-                          MPI_Comm comm)
-{
-  int retval = MPI_SUCCESS;
-  int block_dsize;
-  int rank;
+MPI_CALL_IMPLEM(int, MPI_Cancel, (MPI_Request* request)) {
+   return not_yet_implemented();
+}
 
-  smpi_bench_end();
+MPI_CALL_IMPLEM(int, MPI_Buffer_attach, (void* buffer, int size)) {
+   return not_yet_implemented();
+}
 
-  rank = smpi_mpi_comm_rank(comm);
-  block_dsize = datatype->size * sendcount;
-  INFO2("[%d] optimized alltoall() called. Block size sent to each rank=%d.\n",rank,block_dsize);
+MPI_CALL_IMPLEM(int, MPI_Buffer_detach, (void* buffer, int* size)) {
+   return not_yet_implemented();
+}
 
-  if ((block_dsize < 200) && (comm->size > 12)) {
-           retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
-                                 recvbuf, recvcount, recvtype, comm);
+MPI_CALL_IMPLEM(int, MPI_Testsome, (int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses)) {
+   return not_yet_implemented();
+}
 
-  } else if (block_dsize < 3000) {
-           retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, datatype,
-                                 recvbuf, recvcount, recvtype, comm);
-  } else {
+MPI_CALL_IMPLEM(int, MPI_Comm_test_inter, (MPI_Comm comm, int* flag)) {
+   return not_yet_implemented();
+}
 
-  retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
-                                 recvbuf, recvcount, recvtype, comm);
-  }
+MPI_CALL_IMPLEM(int, MPI_Unpack, (void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm)) {
+   return not_yet_implemented();
+}
 
-  smpi_bench_begin();
+MPI_CALL_IMPLEM(int, MPI_Type_commit, (MPI_Datatype* datatype)) {
+   return not_yet_implemented();
+}
 
-  return retval;
+MPI_CALL_IMPLEM(int, MPI_Type_hindexed, (int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
 }
 
+MPI_CALL_IMPLEM(int, MPI_Type_hvector, (int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
+}
 
+MPI_CALL_IMPLEM(int, MPI_Type_indexed, (int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
+}
 
+MPI_CALL_IMPLEM(int, MPI_Type_struct, (int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
+}
 
-// used by comm_split to sort ranks based on key values
-int smpi_compare_rankkeys(const void *a, const void *b);
-int smpi_compare_rankkeys(const void *a, const void *b)
-{
-  int *x = (int *) a;
-  int *y = (int *) b;
+MPI_CALL_IMPLEM(int, MPI_Type_vector, (int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype)) {
+   return not_yet_implemented();
+}
 
-  if (x[1] < y[1])
-    return -1;
+MPI_CALL_IMPLEM(int, MPI_Ssend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm)) {
+   return not_yet_implemented();
+}
 
-  if (x[1] == y[1]) {
-    if (x[0] < y[0])
-      return -1;
-    if (x[0] == y[0])
-      return 0;
-    return 1;
-  }
+MPI_CALL_IMPLEM(int, MPI_Ssend_init, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
+}
 
-  return 1;
+MPI_CALL_IMPLEM(int, MPI_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();
 }
 
-int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
-                        MPI_Comm * comm_out)
-{
-  int retval = MPI_SUCCESS;
+MPI_CALL_IMPLEM(int, MPI_Intercomm_merge, (MPI_Comm comm, int high, MPI_Comm* comm_out)) {
+   return not_yet_implemented();
+}
 
-  int index, rank;
-  smpi_mpi_request_t request;
-  int colorkey[2];
-  smpi_mpi_status_t status;
+MPI_CALL_IMPLEM(int, MPI_Bsend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm)) {
+   return not_yet_implemented();
+}
 
-  smpi_bench_end();
+MPI_CALL_IMPLEM(int, MPI_Bsend_init, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
+}
 
-  // FIXME: need to test parameters
+MPI_CALL_IMPLEM(int, MPI_Ibsend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
+}
 
-  index = smpi_process_index();
-  rank = comm->index_to_rank_map[index];
+MPI_CALL_IMPLEM(int, MPI_Comm_remote_group, (MPI_Comm comm, MPI_Group* group)) {
+   return not_yet_implemented();
+}
 
-  // default output
-  comm_out = NULL;
+MPI_CALL_IMPLEM(int, MPI_Comm_remote_size, (MPI_Comm comm, int* size)) {
+   return not_yet_implemented();
+}
 
-  // root node does most of the real work
-  if (0 == rank) {
-    int colormap[comm->size];
-    int keymap[comm->size];
-    int rankkeymap[comm->size * 2];
-    int i, j;
-    smpi_mpi_communicator_t tempcomm = NULL;
-    int count;
-    int indextmp;
+MPI_CALL_IMPLEM(int, MPI_Issend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
+}
 
-    colormap[0] = color;
-    keymap[0] = key;
+MPI_CALL_IMPLEM(int, MPI_Probe, (int source, int tag, MPI_Comm comm, MPI_Status* status)) {
+   return not_yet_implemented();
+}
 
-    // FIXME: use scatter/gather or similar instead of individual comms
-    for (i = 1; i < comm->size; i++) {
-      retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
-                                   rank, MPI_ANY_TAG, comm, &request);
-      smpi_mpi_irecv(request);
-      smpi_mpi_wait(request, &status);
-      colormap[status.MPI_SOURCE] = colorkey[0];
-      keymap[status.MPI_SOURCE] = colorkey[1];
-      xbt_mallocator_release(smpi_global->request_mallocator, request);
-    }
+MPI_CALL_IMPLEM(int, MPI_Attr_delete, (MPI_Comm comm, int keyval)) {
+   return not_yet_implemented();
+}
 
-    for (i = 0; i < comm->size; i++) {
-      if (MPI_UNDEFINED == colormap[i]) {
-        continue;
-      }
-      // make a list of nodes with current color and sort by keys
-      count = 0;
-      for (j = i; j < comm->size; j++) {
-        if (colormap[i] == colormap[j]) {
-          colormap[j] = MPI_UNDEFINED;
-          rankkeymap[count * 2] = j;
-          rankkeymap[count * 2 + 1] = keymap[j];
-          count++;
-        }
-      }
-      qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
-
-      // new communicator
-      tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
-      tempcomm->barrier_count = 0;
-      tempcomm->size = count;
-      tempcomm->barrier_mutex = SIMIX_mutex_init();
-      tempcomm->barrier_cond = SIMIX_cond_init();
-      tempcomm->rank_to_index_map = xbt_new(int, count);
-      tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
-      for (j = 0; j < smpi_global->process_count; j++) {
-        tempcomm->index_to_rank_map[j] = -1;
-      }
-      for (j = 0; j < count; j++) {
-        indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
-        tempcomm->rank_to_index_map[j] = indextmp;
-        tempcomm->index_to_rank_map[indextmp] = j;
-      }
-      for (j = 0; j < count; j++) {
-        if (rankkeymap[j * 2]) {
-          retval = smpi_create_request(&j, 1, MPI_INT, 0,
-                                       rankkeymap[j * 2], 0, comm, &request);
-          request->data = tempcomm;
-          smpi_mpi_isend(request);
-          smpi_mpi_wait(request, &status);
-          xbt_mallocator_release(smpi_global->request_mallocator, request);
-        } else {
-          *comm_out = tempcomm;
-        }
-      }
-    }
-  } else {
-    colorkey[0] = color;
-    colorkey[1] = key;
-    retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
-                                 &request);
-    smpi_mpi_isend(request);
-    smpi_mpi_wait(request, &status);
-    xbt_mallocator_release(smpi_global->request_mallocator, request);
-    if (MPI_UNDEFINED != color) {
-      retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
-                                   &request);
-      smpi_mpi_irecv(request);
-      smpi_mpi_wait(request, &status);
-      *comm_out = request->data;
-    }
-  }
+MPI_CALL_IMPLEM(int, MPI_Attr_get, (MPI_Comm comm, int keyval, void* attr_value, int* flag)) {
+   return not_yet_implemented();
+}
 
-  smpi_bench_begin();
+MPI_CALL_IMPLEM(int, MPI_Attr_put, (MPI_Comm comm, int keyval, void* attr_value)) {
+   return not_yet_implemented();
+}
 
-  return retval;
+MPI_CALL_IMPLEM(int, MPI_Rsend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm)) {
+   return not_yet_implemented();
 }
 
-double SMPI_MPI_Wtime(void)
-{
-  return (SIMIX_get_clock());
+MPI_CALL_IMPLEM(int, MPI_Rsend_init, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Irsend, (void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request)) {
+   return not_yet_implemented();
 }
 
+MPI_CALL_IMPLEM(int, MPI_Keyval_create, (MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Keyval_free, (int* keyval)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Test_cancelled, (MPI_Status* status, int* flag)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Pack, (void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Testall, (int count, MPI_Request* requests, int* flag, MPI_Status* statuses)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Get_elements, (MPI_Status* status, MPI_Datatype datatype, int* elements)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Dims_create, (int nnodes, int ndims, int* dims)) {
+   return not_yet_implemented();
+}
+
+MPI_CALL_IMPLEM(int, MPI_Iprobe, (int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status)) {
+   return not_yet_implemented();
+}
 
+MPI_CALL_IMPLEM(int, MPI_Initialized, (int* flag)) {
+   return not_yet_implemented();
+}