Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Fixed return type.
[simgrid.git] / src / smpi / smpi_mpi.c
index ca95953..3bf7d08 100644 (file)
@@ -1,4 +1,8 @@
-/* $Id$tag */
+/* 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"
@@ -12,10 +16,16 @@ XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
 int MPI_Init(int* argc, char*** argv) {
   smpi_process_init(argc, argv);
   smpi_bench_begin(-1, NULL);
+#ifdef HAVE_TRACING
+  TRACE_smpi_init(smpi_process_index());
+#endif
   return MPI_SUCCESS;
 }
 
 int MPI_Finalize(void) {
+#ifdef HAVE_TRACING
+  TRACE_smpi_finalize(smpi_process_index());
+#endif
   smpi_bench_end(-1, NULL);
   smpi_process_destroy();
   return MPI_SUCCESS;
@@ -100,7 +110,7 @@ int MPI_Type_free(MPI_Datatype* datatype) {
   return retval;
 }
 
-int MPI_Type_size(MPI_Datatype datatype, size_t* size) {
+int MPI_Type_size(MPI_Datatype datatype, int* size) {
   int retval;
 
   smpi_bench_end(-1, NULL);
@@ -109,7 +119,7 @@ int MPI_Type_size(MPI_Datatype datatype, size_t* size) {
   } else if(size == NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    *size = smpi_datatype_size(datatype);
+    *size = (int)smpi_datatype_size(datatype);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(-1, NULL);
@@ -681,10 +691,96 @@ int MPI_Comm_free(MPI_Comm* comm) {
   return retval;
 }
 
+int MPI_Send_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
+  int retval;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end(rank, "Send_init");
+  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(rank, "Send_init");
+  return retval;
+}
+
+int MPI_Recv_init(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
+  int retval;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end(rank, "Recv_init");
+  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(rank, "Recv_init");
+  return retval;
+}
+
+int MPI_Start(MPI_Request* request) {
+  int retval;
+  MPI_Comm comm = (*request)->comm;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end(rank, "Start");
+  if(request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_start(*request);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin(rank, "Start");
+  return retval;
+}
+
+int MPI_Startall(int count, MPI_Request* requests) {
+  int retval;
+  MPI_Comm comm = count > 0 && requests ? requests[0]->comm : MPI_COMM_NULL;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end(rank, "Startall");
+  if(requests == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_startall(count, requests);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin(rank, "Startall");
+  return retval;
+}
+
+int MPI_Request_free(MPI_Request* request) {
+  int retval;
+  MPI_Comm comm = (*request)->comm;
+  int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+
+  smpi_bench_end(rank, "Request_free");
+  if(request == NULL) {
+    retval = MPI_ERR_ARG;
+  } else {
+    smpi_mpi_request_free(request);
+    retval = MPI_SUCCESS;
+  }
+  smpi_bench_begin(rank, "Request_free");
+  return retval;
+}
+
 int MPI_Irecv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in (rank, src_traced, rank, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Irecv");
   if(request == NULL) {
     retval = MPI_ERR_ARG;
@@ -695,6 +791,10 @@ int MPI_Irecv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Irecv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out (rank, src_traced, rank, __FUNCTION__);
+  (*request)->recv = 1;
+#endif
   return retval;
 }
 
@@ -702,6 +802,11 @@ int MPI_Isend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  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
   smpi_bench_end(rank, "Isend");
   if(request == NULL) {
     retval = MPI_ERR_ARG;
@@ -712,6 +817,10 @@ int MPI_Isend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Isend");
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out (rank, rank, dst_traced, __FUNCTION__);
+  (*request)->send = 1;
+#endif
   return retval;
 }
 
@@ -719,6 +828,10 @@ int MPI_Recv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
+  TRACE_smpi_ptp_in (rank, src_traced, rank, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Recv");
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -727,6 +840,10 @@ int MPI_Recv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Recv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out (rank, src_traced, rank, __FUNCTION__);
+  TRACE_smpi_recv (rank, src_traced, rank);
+#endif
   return retval;
 }
 
@@ -734,6 +851,11 @@ int MPI_Send(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  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
   smpi_bench_end(rank, "Send");
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -742,6 +864,9 @@ int MPI_Send(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Send");
+#ifdef HAVE_TRACING
+  TRACE_smpi_ptp_out (rank, rank, dst_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -749,6 +874,13 @@ int MPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, i
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  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
   smpi_bench_end(rank, "Sendrecv");
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -759,6 +891,11 @@ int MPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, i
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Sendrecv");
+#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
   return retval;
 }
 
@@ -775,30 +912,6 @@ int MPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, i
   return retval;
 }
 
-int MPI_Get_count(MPI_Status *status, MPI_Datatype datatype, int *count) {
-  int retval;
-/*
- * Returns the number of entries received. (Again, we count entries, each of type datatype, not bytes.)
- * The datatype argument should match the argument provided by the receive call that set the status variable.
- * If the size of the datatype is zero, this routine will return a count of zero.
- * If the amount of data in status is not an exact multiple of the size of datatype
- * (so that count would not be integral), a count of MPI_UNDEFINED is returned instead.
- *
- */
-  smpi_bench_end(-1, NULL); //FIXME
-
-  if( 0==smpi_datatype_size(datatype)) {
-          // also check that the type is 'committed' when we have MPI_Type_commit (s.g. 23/03/21010)
-           retval = MPI_ERR_TYPE;
-  } else {
-           smpi_mpi_get_count(status, datatype, count);
-           retval = MPI_SUCCESS;
-  }
-  smpi_bench_begin(-1, NULL);
-  return retval;
-}
-
-
 int MPI_Test(MPI_Request* request, int* flag, MPI_Status* status) {
   int retval;
   int rank = request && (*request)->comm != MPI_COMM_NULL
@@ -838,6 +951,13 @@ int MPI_Wait(MPI_Request* request, MPI_Status* status) {
              ? smpi_comm_rank((*request)->comm)
              : -1;
 
+#ifdef HAVE_TRACING
+  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
   smpi_bench_end(rank, "Wait");
   if(request == NULL) {
     retval = MPI_ERR_ARG;
@@ -848,12 +968,53 @@ int MPI_Wait(MPI_Request* request, MPI_Status* status) {
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Wait");
+#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
   return retval;
 }
 
 int MPI_Waitany(int count, MPI_Request requests[], int* index, MPI_Status* status) {
   int retval;
 
+#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);
+    }
+  }
+
+  //search for a suitable request to give the rank of current mpi proc
+  MPI_Request req = NULL;
+  for (i = 0; i < count && req == NULL; i++) {
+    req = requests[i];
+  }
+  MPI_Comm comm = (req)->comm;
+  int rank_traced = smpi_comm_rank(comm);
+  TRACE_smpi_ptp_in (rank_traced, -1, -1, __FUNCTION__);
+#endif
   smpi_bench_end(-1, NULL); //FIXME
   if(index == NULL) {
     retval = MPI_ERR_ARG;
@@ -862,13 +1023,72 @@ int MPI_Waitany(int count, MPI_Request requests[], int* index, MPI_Status* statu
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(-1, NULL);
+#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
   return retval;
 }
 
 int MPI_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), 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);
+  }
+
+//  find my rank inside one of MPI_Comm's of the requests
+  MPI_Request req = NULL;
+  for (i = 0; i < count && req == NULL; i++) {
+    req = requests[i];
+  }
+  MPI_Comm comm = (req)->comm;
+  int rank_traced = smpi_comm_rank(comm);
+  TRACE_smpi_ptp_in (rank_traced, -1, -1, __FUNCTION__);
+#endif
   smpi_bench_end(-1, NULL); //FIXME
   smpi_mpi_waitall(count, requests, status);
   smpi_bench_begin(-1, NULL);
+#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
   return MPI_SUCCESS;
 }
 
@@ -890,6 +1110,10 @@ int MPI_Bcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm co
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Bcast");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -898,6 +1122,9 @@ int MPI_Bcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm co
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Bcast");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -905,6 +1132,9 @@ int MPI_Barrier(MPI_Comm comm) {
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Barrier");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -913,6 +1143,9 @@ int MPI_Barrier(MPI_Comm comm) {
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Barrier");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -920,6 +1153,10 @@ int MPI_Gather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbu
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Gather");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -930,6 +1167,9 @@ int MPI_Gather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbu
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Gather");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -937,6 +1177,10 @@ int MPI_Gatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvb
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Gatherv");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -949,6 +1193,9 @@ int MPI_Gatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvb
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Gatherv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -956,6 +1203,9 @@ int MPI_Allgather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* rec
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Allgather");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -966,6 +1216,9 @@ int MPI_Allgather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* rec
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Allgather");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -973,6 +1226,9 @@ int MPI_Allgatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* re
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Allgatherv");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -985,6 +1241,9 @@ int MPI_Allgatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* re
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Allgatherv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -992,6 +1251,10 @@ int MPI_Scatter(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvb
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Scatter");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1002,6 +1265,9 @@ int MPI_Scatter(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvb
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Scatter");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1009,6 +1275,10 @@ int MPI_Scatterv(void* sendbuf, int* sendcounts, int* displs, MPI_Datatype sendt
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Scatterv");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1021,6 +1291,9 @@ int MPI_Scatterv(void* sendbuf, int* sendcounts, int* displs, MPI_Datatype sendt
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Scatterv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1028,6 +1301,10 @@ int MPI_Reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, M
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
+  TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Reduce");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1038,6 +1315,9 @@ int MPI_Reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, M
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Reduce");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1045,6 +1325,9 @@ int MPI_Allreduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Allreduce");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1057,6 +1340,9 @@ int MPI_Allreduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Allreduce");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1064,6 +1350,9 @@ int MPI_Scan(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Scan");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1076,6 +1365,9 @@ int MPI_Scan(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Scan");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1084,6 +1376,9 @@ int MPI_Reduce_scatter(void* sendbuf, void* recvbuf, int* recvcounts, MPI_Dataty
   int* displs;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Reduce_scatter");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1109,6 +1404,9 @@ int MPI_Reduce_scatter(void* sendbuf, void* recvbuf, int* recvcounts, MPI_Dataty
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin(rank, "Reduce_scatter");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1116,6 +1414,9 @@ int MPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recv
   int retval, size, sendsize;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Alltoall");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1133,6 +1434,9 @@ int MPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recv
     }
   }
   smpi_bench_begin(rank, "Alltoall");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
@@ -1140,6 +1444,9 @@ int MPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype s
   int retval;
   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
 
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
+#endif
   smpi_bench_end(rank, "Alltoallv");
   if(comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
@@ -1151,12 +1458,16 @@ int MPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype s
     retval = smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm);
   }
   smpi_bench_begin(rank, "Alltoallv");
+#ifdef HAVE_TRACING
+  TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
+#endif
   return retval;
 }
 
 
-int MPI_Get_processor_name( char *name, int *resultlen ) {
+int MPI_Get_processor_name(char* name, int* resultlen) {
   int retval = MPI_SUCCESS;
+
   smpi_bench_end(-1, NULL);
   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);
@@ -1165,3 +1476,25 @@ int MPI_Get_processor_name( char *name, int *resultlen ) {
   return retval;
 }
 
+int MPI_Get_count(MPI_Status* status, MPI_Datatype datatype, int* count) {
+  int retval = MPI_SUCCESS;
+  size_t size;
+
+  smpi_bench_end(-1, NULL);
+  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(-1, NULL);
+  return retval;
+}