Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[TRACING] Rename TIT action waitAll -> waitall
[simgrid.git] / src / smpi / bindings / smpi_pmpi_request.cpp
index 88a59c8..19e2f99 100644 (file)
@@ -15,11 +15,10 @@ static int getPid(MPI_Comm, int);
 static int getPid(MPI_Comm comm, int id)
 {
   simgrid::s4u::ActorPtr actor = comm->group()->actor(id);
-  return (actor == nullptr) ? MPI_UNDEFINED : actor->getPid();
+  return (actor == nullptr) ? MPI_UNDEFINED : actor->get_pid();
 }
 
 /* PMPI User level calls */
-extern "C" { // Obviously, the C MPI interface should use the C linkage
 
 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
 {
@@ -90,6 +89,11 @@ int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int ta
   return retval;
 }
 
+/*
+ * This function starts a request returned by init functions such as
+ * MPI_Send_init(), MPI_Ssend_init (see above), and friends.
+ * They should already have performed sanity checks.
+ */
 int PMPI_Start(MPI_Request * request)
 {
   int retval = 0;
@@ -98,8 +102,22 @@ int PMPI_Start(MPI_Request * request)
   if (request == nullptr || *request == MPI_REQUEST_NULL) {
     retval = MPI_ERR_REQUEST;
   } else {
-    (*request)->start();
+    MPI_Request req = *request;
+    int my_proc_id = (req->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
+    TRACE_smpi_comm_in(my_proc_id, __func__,
+                       new simgrid::instr::Pt2PtTIData("Start", req->dst(),
+                                                       req->size(),
+                                                       req->tag(), 
+                                                       simgrid::smpi::Datatype::encode(req->type())));
+    if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_SEND)
+      TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
+
+    req->start();
+
+    if (not TRACE_smpi_view_internals() && req->flags() & MPI_REQ_RECV)
+      TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
     retval = MPI_SUCCESS;
+    TRACE_smpi_comm_out(my_proc_id);
   }
   smpi_bench_begin();
   return retval;
@@ -119,7 +137,25 @@ int PMPI_Startall(int count, MPI_Request * requests)
       }
     }
     if(retval != MPI_ERR_REQUEST) {
+      int my_proc_id = simgrid::s4u::this_actor::get_pid();
+      TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("Startall"));
+      MPI_Request req = MPI_REQUEST_NULL;
+      if (not TRACE_smpi_view_internals())
+        for (int i = 0; i < count; i++) {
+          req = requests[i];
+          if (req->flags() & MPI_REQ_SEND)
+            TRACE_smpi_send(my_proc_id, my_proc_id, getPid(req->comm(), req->dst()), req->tag(), req->size());
+        }
+
       simgrid::smpi::Request::startall(count, requests);
+
+      if (not TRACE_smpi_view_internals())
+        for (int i = 0; i < count; i++) {
+          req = requests[i];
+          if (req->flags() & MPI_REQ_RECV)
+            TRACE_smpi_recv(getPid(req->comm(), req->src()), my_proc_id, req->tag());
+        }
+      TRACE_smpi_comm_out(my_proc_id);
     }
   }
   smpi_bench_begin();
@@ -164,12 +200,12 @@ int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MP
     retval = MPI_ERR_TAG;
   } else {
 
-    int my_proc_id = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("Irecv", src,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
 
     *request = simgrid::smpi::Request::irecv(buf, count, datatype, src, tag, comm);
     retval = MPI_SUCCESS;
@@ -205,12 +241,12 @@ int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MP
   } else if(tag<0 && tag !=  MPI_ANY_TAG){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
     int trace_dst = getPid(comm, dst);
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("Isend", dst,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
 
     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
 
@@ -247,12 +283,12 @@ int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, M
   } else if(tag<0 && tag !=  MPI_ANY_TAG){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
     int trace_dst = getPid(comm, dst);
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("ISsend", dst,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
     TRACE_smpi_send(my_proc_id, my_proc_id, trace_dst, tag, count * datatype->size());
 
     *request = simgrid::smpi::Request::issend(buf, count, datatype, dst, tag, comm);
@@ -275,8 +311,10 @@ int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else if (src == MPI_PROC_NULL) {
-    simgrid::smpi::Status::empty(status);
-    status->MPI_SOURCE = MPI_PROC_NULL;
+    if(status != MPI_STATUS_IGNORE){
+      simgrid::smpi::Status::empty(status);
+      status->MPI_SOURCE = MPI_PROC_NULL;
+    }
     retval = MPI_SUCCESS;
   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
     retval = MPI_ERR_RANK;
@@ -287,22 +325,25 @@ int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI
   } else if(tag<0 && tag !=  MPI_ANY_TAG){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id         = simgrid::s4u::Actor::self()->getPid();
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("recv", src,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
 
     simgrid::smpi::Request::recv(buf, count, datatype, src, tag, comm, status);
     retval = MPI_SUCCESS;
 
     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
-    if (status != MPI_STATUS_IGNORE) {
-      int src_traced = getPid(comm, status->MPI_SOURCE);
-      if (not TRACE_smpi_view_internals()) {
-        TRACE_smpi_recv(src_traced, my_proc_id, tag);
-      }
+    int src_traced=0;
+    if (status != MPI_STATUS_IGNORE) 
+      src_traced = getPid(comm, status->MPI_SOURCE);
+    else
+      src_traced = getPid(comm, src);
+    if (not TRACE_smpi_view_internals()) {
+      TRACE_smpi_recv(src_traced, my_proc_id, tag);
     }
+    
     TRACE_smpi_comm_out(my_proc_id);
   }
 
@@ -329,12 +370,12 @@ int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI
   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id         = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id         = simgrid::s4u::this_actor::get_pid();
     int dst_traced         = getPid(comm, dst);
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("send", dst,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
     if (not TRACE_smpi_view_internals()) {
       TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
     }
@@ -367,12 +408,12 @@ int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MP
   } else if(tag<0 && tag !=  MPI_ANY_TAG){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id         = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id         = simgrid::s4u::this_actor::get_pid();
     int dst_traced         = getPid(comm, dst);
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::Pt2PtTIData("Ssend", dst,
                                                        datatype->is_replayable() ? count : count * datatype->size(),
-                                                       encode_datatype(datatype)));
+                                                       tag, simgrid::smpi::Datatype::encode(datatype)));
     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, tag, count * datatype->size());
 
     simgrid::smpi::Request::ssend(buf, count, datatype, dst, tag, comm);
@@ -396,20 +437,17 @@ int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst,
     retval = MPI_ERR_COMM;
   } else if (not sendtype->is_valid() || not recvtype->is_valid()) {
     retval = MPI_ERR_TYPE;
-  } else if (src == MPI_PROC_NULL && dst == MPI_PROC_NULL) {
-    if(status!=MPI_STATUS_IGNORE){
-      simgrid::smpi::Status::empty(status);
-      status->MPI_SOURCE = MPI_PROC_NULL;
-    }
-    retval             = MPI_SUCCESS;
-  }else if (src == MPI_PROC_NULL){
+  } else if (src == MPI_PROC_NULL) {
     if(status!=MPI_STATUS_IGNORE){
       simgrid::smpi::Status::empty(status);
       status->MPI_SOURCE = MPI_PROC_NULL;
     }
-    return PMPI_Send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
+    if(dst != MPI_PROC_NULL)
+      simgrid::smpi::Request::send(sendbuf, sendcount, sendtype, dst, sendtag, comm);
+    retval = MPI_SUCCESS;
   }else if (dst == MPI_PROC_NULL){
-    return PMPI_Recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
+    simgrid::smpi::Request::recv(recvbuf, recvcount, recvtype, src, recvtag, comm, status);
+    retval = MPI_SUCCESS;
   }else if (dst >= comm->group()->size() || dst <0 ||
       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
     retval = MPI_ERR_RANK;
@@ -419,7 +457,7 @@ int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst,
   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
     retval = MPI_ERR_TAG;
   } else {
-    int my_proc_id         = simgrid::s4u::Actor::self()->getPid();
+    int my_proc_id         = simgrid::s4u::this_actor::get_pid();
     int dst_traced         = getPid(comm, dst);
     int src_traced         = getPid(comm, src);
 
@@ -428,11 +466,11 @@ int PMPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst,
     std::vector<int>* src_hack = new std::vector<int>;
     dst_hack->push_back(dst_traced);
     src_hack->push_back(src_traced);
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__,
+    TRACE_smpi_comm_in(my_proc_id, __func__,
                        new simgrid::instr::VarCollTIData(
                            "sendRecv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
                            dst_hack, recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), src_hack,
-                           encode_datatype(sendtype), encode_datatype(recvtype)));
+                           simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
 
     TRACE_smpi_send(my_proc_id, my_proc_id, dst_traced, sendtag, sendcount * sendtype->size());
 
@@ -476,17 +514,19 @@ int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
   if (request == nullptr || flag == nullptr) {
     retval = MPI_ERR_ARG;
   } else if (*request == MPI_REQUEST_NULL) {
-    *flag= true;
-    simgrid::smpi::Status::empty(status);
+    if (status != MPI_STATUS_IGNORE){
+      *flag= true;
+      simgrid::smpi::Status::empty(status);
+    }
     retval = MPI_SUCCESS;
   } else {
-    int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::Actor::self()->getPid() : -1;
-
-    TRACE_smpi_testing_in(my_proc_id);
+    int my_proc_id = ((*request)->comm() != MPI_COMM_NULL) ? simgrid::s4u::this_actor::get_pid() : -1;
 
+    TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("test"));
+    
     *flag = simgrid::smpi::Request::test(request,status);
 
-    TRACE_smpi_testing_out(my_proc_id);
+    TRACE_smpi_comm_out(my_proc_id);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -501,7 +541,10 @@ int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_S
   if (index == nullptr || flag == nullptr) {
     retval = MPI_ERR_ARG;
   } else {
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testany"));
     *flag = simgrid::smpi::Request::testany(count, requests, index, status);
+    TRACE_smpi_comm_out(my_proc_id);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
@@ -516,24 +559,45 @@ int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* status
   if (flag == nullptr) {
     retval = MPI_ERR_ARG;
   } else {
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testall"));
     *flag = simgrid::smpi::Request::testall(count, requests, statuses);
+    TRACE_smpi_comm_out(my_proc_id);
     retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
   return retval;
 }
 
+int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
+{
+  int retval = 0;
+
+  smpi_bench_end();
+  if (outcount == nullptr) {
+    retval = MPI_ERR_ARG;
+  } else {
+    int my_proc_id = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(my_proc_id, __func__, new simgrid::instr::NoOpTIData("testsome"));
+    *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
+    TRACE_smpi_comm_out(my_proc_id);
+    retval    = MPI_SUCCESS;
+  }
+  smpi_bench_begin();
+  return retval;
+}
+
 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
   int retval = 0;
   smpi_bench_end();
 
-  if (status == nullptr) {
-    retval = MPI_ERR_ARG;
-  } else if (comm == MPI_COMM_NULL) {
+  if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else if (source == MPI_PROC_NULL) {
-    simgrid::smpi::Status::empty(status);
-    status->MPI_SOURCE = MPI_PROC_NULL;
+    if (status != MPI_STATUS_IGNORE){
+      simgrid::smpi::Status::empty(status);
+      status->MPI_SOURCE = MPI_PROC_NULL;
+    }
     retval = MPI_SUCCESS;
   } else {
     simgrid::smpi::Request::probe(source, tag, comm, status);
@@ -553,8 +617,10 @@ int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* statu
     retval = MPI_ERR_COMM;
   } else if (source == MPI_PROC_NULL) {
     *flag=true;
-    simgrid::smpi::Status::empty(status);
-    status->MPI_SOURCE = MPI_PROC_NULL;
+    if (status != MPI_STATUS_IGNORE){
+      simgrid::smpi::Status::empty(status);
+      status->MPI_SOURCE = MPI_PROC_NULL;
+    }
     retval = MPI_SUCCESS;
   } else {
     simgrid::smpi::Request::iprobe(source, tag, comm, flag, status);
@@ -573,9 +639,9 @@ static void trace_smpi_recv_helper(MPI_Request* request, MPI_Status* status)
     int src_traced = req->src();
     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
     int dst_traced = req->dst();
-    if (req->flags() & RECV) { // Is this request a wait for RECV?
+    if (req->flags() & MPI_REQ_RECV) { // Is this request a wait for RECV?
       if (src_traced == MPI_ANY_SOURCE)
-        src_traced = (status != MPI_STATUSES_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
+        src_traced = (status != MPI_STATUS_IGNORE) ? req->comm()->group()->rank(status->MPI_SOURCE) : req->src();
       TRACE_smpi_recv(src_traced, dst_traced, req->tag());
     }
   }
@@ -594,18 +660,27 @@ int PMPI_Wait(MPI_Request * request, MPI_Status * status)
   } else if (*request == MPI_REQUEST_NULL) {
     retval = MPI_SUCCESS;
   } else {
+    //for tracing, save the handle which might get overriden before we can use the helper on it
+    MPI_Request savedreq = *request;
+    if (savedreq != MPI_REQUEST_NULL && not(savedreq->flags() & MPI_REQ_FINISHED))
+      savedreq->ref();//don't erase te handle in Request::wait, we'll need it later
+    else
+      savedreq = MPI_REQUEST_NULL;
+
     int my_proc_id = (*request)->comm() != MPI_COMM_NULL
-                         ? simgrid::s4u::Actor::self()->getPid()
+                         ? simgrid::s4u::this_actor::get_pid()
                          : -1; // TODO: cheinrich: Check if this correct or if it should be MPI_UNDEFINED
-
-    TRACE_smpi_comm_in(my_proc_id, __FUNCTION__, new simgrid::instr::NoOpTIData("wait"));
+    TRACE_smpi_comm_in(my_proc_id, __func__,
+                       new simgrid::instr::WaitTIData((*request)->src(), (*request)->dst(), (*request)->tag()));
 
     simgrid::smpi::Request::wait(request, status);
     retval = MPI_SUCCESS;
 
     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
     TRACE_smpi_comm_out(my_proc_id);
-    trace_smpi_recv_helper(request, status);
+    trace_smpi_recv_helper(&savedreq, status);
+    if (savedreq != MPI_REQUEST_NULL)
+      simgrid::smpi::Request::unref(&savedreq);
   }
 
   smpi_bench_begin();
@@ -621,17 +696,29 @@ int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
     return MPI_SUCCESS;
 
   smpi_bench_end();
+  //for tracing, save the handles which might get overriden before we can use the helper on it
+  std::vector<MPI_Request> savedreqs(requests, requests + count);
+  for (MPI_Request& req : savedreqs) {
+    if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
+      req->ref();
+    else
+      req = MPI_REQUEST_NULL;
+  }
 
-  int rank_traced = simgrid::s4u::Actor::self()->getPid(); // FIXME: In PMPI_Wait, we check if the comm is null?
-  TRACE_smpi_comm_in(rank_traced, __FUNCTION__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
+  int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
+  TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitAny", static_cast<double>(count)));
 
   *index = simgrid::smpi::Request::waitany(count, requests, status);
 
   if(*index!=MPI_UNDEFINED){
-    trace_smpi_recv_helper(&requests[*index], status);
+    trace_smpi_recv_helper(&savedreqs[*index], status);
     TRACE_smpi_comm_out(rank_traced);
   }
 
+  for (MPI_Request& req : savedreqs)
+    if (req != MPI_REQUEST_NULL)
+      simgrid::smpi::Request::unref(&req);
+
   smpi_bench_begin();
   return MPI_SUCCESS;
 }
@@ -640,16 +727,29 @@ int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
 {
   smpi_bench_end();
 
-  int rank_traced = simgrid::s4u::Actor::self()->getPid(); // FIXME: In PMPI_Wait, we check if the comm is null?
-  TRACE_smpi_comm_in(rank_traced, __FUNCTION__, new simgrid::instr::CpuTIData("waitAll", static_cast<double>(count)));
+  //for tracing, save the handles which might get overriden before we can use the helper on it
+  std::vector<MPI_Request> savedreqs(requests, requests + count);
+  for (MPI_Request& req : savedreqs) {
+    if (req != MPI_REQUEST_NULL && not(req->flags() & MPI_REQ_FINISHED))
+      req->ref();
+    else
+      req = MPI_REQUEST_NULL;
+  }
+
+  int rank_traced = simgrid::s4u::this_actor::get_pid(); // FIXME: In PMPI_Wait, we check if the comm is null?
+  TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("waitall", static_cast<double>(count)));
 
   int retval = simgrid::smpi::Request::waitall(count, requests, status);
 
   for (int i = 0; i < count; i++) {
-    trace_smpi_recv_helper(&requests[i], &status[i]);
+    trace_smpi_recv_helper(&savedreqs[i], status!=MPI_STATUSES_IGNORE ? &status[i]: MPI_STATUS_IGNORE);
   }
   TRACE_smpi_comm_out(rank_traced);
 
+  for (MPI_Request& req : savedreqs)
+    if (req != MPI_REQUEST_NULL)
+      simgrid::smpi::Request::unref(&req);
+
   smpi_bench_begin();
   return retval;
 }
@@ -669,21 +769,30 @@ int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indic
   return retval;
 }
 
-int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
+int PMPI_Cancel(MPI_Request* request)
 {
   int retval = 0;
 
   smpi_bench_end();
-  if (outcount == nullptr) {
-    retval = MPI_ERR_ARG;
+  if (*request == MPI_REQUEST_NULL) {
+    retval = MPI_ERR_REQUEST;
   } else {
-    *outcount = simgrid::smpi::Request::testsome(incount, requests, indices, status);
-    retval    = MPI_SUCCESS;
+    (*request)->cancel();
+    retval = MPI_SUCCESS;
   }
   smpi_bench_begin();
   return retval;
 }
 
+int PMPI_Test_cancelled(MPI_Status* status, int* flag){
+  if(status==MPI_STATUS_IGNORE){
+    *flag=0;
+    return MPI_ERR_ARG;
+  }
+  *flag=simgrid::smpi::Status::cancelled(status);
+  return MPI_SUCCESS;  
+}
+
 MPI_Request PMPI_Request_f2c(MPI_Fint request){
   return static_cast<MPI_Request>(simgrid::smpi::Request::f2c(request));
 }
@@ -691,5 +800,3 @@ MPI_Request PMPI_Request_f2c(MPI_Fint request){
 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
   return request->c2f();
 }
-
-}