Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[TRACING] Rename TIT action reduceScatter -> reducescatter
[simgrid.git] / src / smpi / bindings / smpi_pmpi_coll.cpp
index 8904505..6fa2231 100644 (file)
@@ -1,4 +1,4 @@
-/* Copyright (c) 2007-2017. The SimGrid Team. All rights reserved.          */
+/* Copyright (c) 2007-2018. 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. */
@@ -13,7 +13,6 @@
 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
 
 /* PMPI User level calls */
-extern "C" { // Obviously, the C MPI interface should use the C linkage
 
 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
 {
@@ -26,11 +25,11 @@ int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm c
   } else if (not datatype->is_valid()) {
     retval = MPI_ERR_ARG;
   } else {
-    int rank = simgrid::s4u::Actor::self()->getPid();
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
+    int rank = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(rank, __func__,
                        new simgrid::instr::CollTIData("bcast", root, -1.0,
                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
-                                                      encode_datatype(datatype), ""));
+                                                      simgrid::smpi::Datatype::encode(datatype), ""));
     if (comm->size() > 1)
       simgrid::smpi::Colls::bcast(buf, count, datatype, root, comm);
     retval = MPI_SUCCESS;
@@ -50,8 +49,8 @@ int PMPI_Barrier(MPI_Comm comm)
   if (comm == MPI_COMM_NULL) {
     retval = MPI_ERR_COMM;
   } else {
-    int rank = simgrid::s4u::Actor::self()->getPid();
-    TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::NoOpTIData("barrier"));
+    int rank = simgrid::s4u::this_actor::get_pid();
+    TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::NoOpTIData("barrier"));
 
     simgrid::smpi::Colls::barrier(comm);
 
@@ -90,14 +89,14 @@ int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbu
       sendtmpcount=0;
       sendtmptype=recvtype;
     }
-    int rank = simgrid::s4u::Actor::self()->getPid();
+    int rank = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::CollTIData(
-                           "gather", root, -1.0,
-                           sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
-                           (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
-                           encode_datatype(sendtmptype), encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(
+        rank, __func__,
+        new simgrid::instr::CollTIData(
+            "gather", root, -1.0, sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
+            (comm->rank() != root || recvtype->is_replayable()) ? recvcount : recvcount * recvtype->size(),
+            simgrid::smpi::Datatype::encode(sendtmptype), simgrid::smpi::Datatype::encode(recvtype)));
 
     simgrid::smpi::Colls::gather(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
 
@@ -123,7 +122,7 @@ int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recv
     retval = MPI_ERR_TYPE;
   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
     retval = MPI_ERR_COUNT;
-  } else if (recvcounts == nullptr || displs == nullptr) {
+  } else if ((comm->rank() == root) && (recvcounts == nullptr || displs == nullptr)) {
     retval = MPI_ERR_ARG;
   } else {
     char* sendtmpbuf = static_cast<char*>(sendbuf);
@@ -134,7 +133,7 @@ int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recv
       sendtmptype=recvtype;
     }
 
-    int rank         = simgrid::s4u::Actor::self()->getPid();
+    int rank         = simgrid::s4u::this_actor::get_pid();
     int dt_size_recv = recvtype->is_replayable() ? 1 : recvtype->size();
 
     std::vector<int>* trace_recvcounts = new std::vector<int>;
@@ -143,11 +142,12 @@ int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recv
         trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
+    TRACE_smpi_comm_in(rank, __func__,
                        new simgrid::instr::VarCollTIData(
                            "gatherV", root,
                            sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(), nullptr,
-                           dt_size_recv, trace_recvcounts, encode_datatype(sendtmptype), encode_datatype(recvtype)));
+                           dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtmptype),
+                           simgrid::smpi::Datatype::encode(recvtype)));
 
     retval = simgrid::smpi::Colls::gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
     TRACE_smpi_comm_out(rank);
@@ -178,13 +178,13 @@ int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
       sendcount=recvcount;
       sendtype=recvtype;
     }
-    int rank = simgrid::s4u::Actor::self()->getPid();
+    int rank = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::CollTIData("allGather", -1, -1.0,
-                                                      sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
-                                                      recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
-                                                      encode_datatype(sendtype), encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::CollTIData(
+                           "allGather", -1, -1.0, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
+                           recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
+                           simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
 
     simgrid::smpi::Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
     TRACE_smpi_comm_out(rank);
@@ -215,17 +215,18 @@ int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
       sendcount=recvcounts[comm->rank()];
       sendtype=recvtype;
     }
-    int rank               = simgrid::s4u::Actor::self()->getPid();
+    int rank               = simgrid::s4u::this_actor::get_pid();
     int dt_size_recv       = recvtype->is_replayable() ? 1 : recvtype->size();
 
     std::vector<int>* trace_recvcounts = new std::vector<int>;
     for (int i = 0; i < comm->size(); i++) // copy data to avoid bad free
       trace_recvcounts->push_back(recvcounts[i] * dt_size_recv);
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
+    TRACE_smpi_comm_in(rank, __func__,
                        new simgrid::instr::VarCollTIData(
-                           "allGatherV", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(), nullptr,
-                           dt_size_recv, trace_recvcounts, encode_datatype(sendtype), encode_datatype(recvtype)));
+                           "allgatherv", -1, sendtype->is_replayable() ? sendcount : sendcount * sendtype->size(),
+                           nullptr, dt_size_recv, trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
+                           simgrid::smpi::Datatype::encode(recvtype)));
 
     simgrid::smpi::Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
     retval = MPI_SUCCESS;
@@ -257,14 +258,15 @@ int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
       recvtype  = sendtype;
       recvcount = sendcount;
     }
-    int rank = simgrid::s4u::Actor::self()->getPid();
+    int rank = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::CollTIData(
-                           "scatter", root, -1.0,
-                           (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
-                           recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), encode_datatype(sendtype),
-                           encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(
+        rank, __func__,
+        new simgrid::instr::CollTIData(
+            "scatter", root, -1.0,
+            (comm->rank() != root || sendtype->is_replayable()) ? sendcount : sendcount * sendtype->size(),
+            recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
+            simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
 
     simgrid::smpi::Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
     retval = MPI_SUCCESS;
@@ -294,7 +296,7 @@ int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
       recvtype  = sendtype;
       recvcount = sendcounts[comm->rank()];
     }
-    int rank               = simgrid::s4u::Actor::self()->getPid();
+    int rank               = simgrid::s4u::this_actor::get_pid();
     int dt_size_send       = sendtype->is_replayable() ? 1 : sendtype->size();
 
     std::vector<int>* trace_sendcounts = new std::vector<int>;
@@ -303,10 +305,11 @@ int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
         trace_sendcounts->push_back(sendcounts[i] * dt_size_send);
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::VarCollTIData(
-                                               "scatterV", root, dt_size_send, trace_sendcounts,
-                                               recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), nullptr,
-                                               encode_datatype(sendtype), encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::VarCollTIData(
+                           "scatterv", root, dt_size_send, trace_sendcounts,
+                           recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(), nullptr,
+                           simgrid::smpi::Datatype::encode(sendtype), simgrid::smpi::Datatype::encode(recvtype)));
 
     retval = simgrid::smpi::Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
 
@@ -328,12 +331,12 @@ int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
   } else if (not datatype->is_valid() || op == MPI_OP_NULL) {
     retval = MPI_ERR_ARG;
   } else {
-    int rank = simgrid::s4u::Actor::self()->getPid();
+    int rank = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
+    TRACE_smpi_comm_in(rank, __func__,
                        new simgrid::instr::CollTIData("reduce", root, 0,
                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
-                                                      encode_datatype(datatype), ""));
+                                                      simgrid::smpi::Datatype::encode(datatype), ""));
 
     simgrid::smpi::Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
 
@@ -378,12 +381,12 @@ int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatyp
       sendtmpbuf = static_cast<char*>(xbt_malloc(count*datatype->get_extent()));
       simgrid::smpi::Datatype::copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
     }
-    int rank = simgrid::s4u::Actor::self()->getPid();
+    int rank = simgrid::s4u::this_actor::get_pid();
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::CollTIData("allReduce", -1, 0,
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::CollTIData("allreduce", -1, 0,
                                                       datatype->is_replayable() ? count : count * datatype->size(), -1,
-                                                      encode_datatype(datatype), ""));
+                                                      simgrid::smpi::Datatype::encode(datatype), ""));
 
     simgrid::smpi::Colls::allreduce(sendtmpbuf, recvbuf, count, datatype, op, comm);
 
@@ -411,15 +414,21 @@ int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MP
   } else if (op == MPI_OP_NULL) {
     retval = MPI_ERR_OP;
   } else {
-    int rank = simgrid::s4u::Actor::self()->getPid();
-
-    TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::Pt2PtTIData(
-                                               "scan", -1, datatype->is_replayable() ? count : count * datatype->size(),
-                                               encode_datatype(datatype)));
+    int rank = simgrid::s4u::this_actor::get_pid();
+    void* sendtmpbuf = sendbuf;
+    if (sendbuf == MPI_IN_PLACE) {
+      sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
+      memcpy(sendtmpbuf, recvbuf, count * datatype->size());
+    }
+    TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
+                                           "scan", -1, datatype->is_replayable() ? count : count * datatype->size(),
+                                           simgrid::smpi::Datatype::encode(datatype)));
 
-    retval = simgrid::smpi::Colls::scan(sendbuf, recvbuf, count, datatype, op, comm);
+    retval = simgrid::smpi::Colls::scan(sendtmpbuf, recvbuf, count, datatype, op, comm);
 
     TRACE_smpi_comm_out(rank);
+    if (sendbuf == MPI_IN_PLACE)
+      xbt_free(sendtmpbuf);
   }
 
   smpi_bench_begin();
@@ -438,16 +447,16 @@ int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
   } else if (op == MPI_OP_NULL) {
     retval = MPI_ERR_OP;
   } else {
-    int rank         = simgrid::s4u::Actor::self()->getPid();
+    int rank         = simgrid::s4u::this_actor::get_pid();
     void* sendtmpbuf = sendbuf;
     if (sendbuf == MPI_IN_PLACE) {
       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::Pt2PtTIData(
-                                               "exscan", -1, datatype->is_replayable() ? count : count * datatype->size(),
-                                               encode_datatype(datatype)));
+    TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::Pt2PtTIData(
+                                           "exscan", -1, datatype->is_replayable() ? count : count * datatype->size(),
+                                           simgrid::smpi::Datatype::encode(datatype)));
 
     retval = simgrid::smpi::Colls::exscan(sendtmpbuf, recvbuf, count, datatype, op, comm);
 
@@ -474,7 +483,7 @@ int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datat
   } else if (recvcounts == nullptr) {
     retval = MPI_ERR_ARG;
   } else {
-    int rank                           = simgrid::s4u::Actor::self()->getPid();
+    int rank                           = simgrid::s4u::this_actor::get_pid();
     std::vector<int>* trace_recvcounts = new std::vector<int>;
     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
     int totalcount    = 0;
@@ -490,9 +499,9 @@ int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datat
       memcpy(sendtmpbuf, recvbuf, totalcount * datatype->size());
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::VarCollTIData("reduceScatter", -1, dt_send_size, nullptr, -1,
-                                                         trace_recvcounts, encode_datatype(datatype), ""));
+    TRACE_smpi_comm_in(rank, __func__, new simgrid::instr::VarCollTIData(
+                                           "reducescatter", -1, dt_send_size, nullptr, -1, trace_recvcounts,
+                                           simgrid::smpi::Datatype::encode(datatype), ""));
 
     simgrid::smpi::Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
     retval = MPI_SUCCESS;
@@ -523,7 +532,7 @@ int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
   } else {
     int count = comm->size();
 
-    int rank                           = simgrid::s4u::Actor::self()->getPid();
+    int rank                           = simgrid::s4u::this_actor::get_pid();
     int dt_send_size                   = datatype->is_replayable() ? 1 : datatype->size();
     std::vector<int>* trace_recvcounts = new std::vector<int>(recvcount * dt_send_size); // copy data to avoid bad free
 
@@ -533,9 +542,9 @@ int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
       memcpy(sendtmpbuf, recvbuf, recvcount * count * datatype->size());
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__,
-                       new simgrid::instr::VarCollTIData("reduceScatter", -1, 0, nullptr, -1, trace_recvcounts,
-                                                         encode_datatype(datatype), ""));
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::VarCollTIData("reducescatter", -1, 0, nullptr, -1, trace_recvcounts,
+                                                         simgrid::smpi::Datatype::encode(datatype), ""));
 
     int* recvcounts = new int[count];
     for (int i      = 0; i < count; i++)
@@ -565,7 +574,7 @@ int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* rec
   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL) {
     retval = MPI_ERR_TYPE;
   } else {
-    int rank                 = simgrid::s4u::Actor::self()->getPid();
+    int rank                 = simgrid::s4u::this_actor::get_pid();
     void* sendtmpbuf         = static_cast<char*>(sendbuf);
     int sendtmpcount         = sendcount;
     MPI_Datatype sendtmptype = sendtype;
@@ -576,12 +585,12 @@ int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* rec
       sendtmptype  = recvtype;
     }
 
-    TRACE_smpi_comm_in(
-        rank, __FUNCTION__,
-        new simgrid::instr::CollTIData("allToAll", -1, -1.0,
-                                       sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
-                                       recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
-                                       encode_datatype(sendtmptype), encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::CollTIData(
+                           "allToAll", -1, -1.0,
+                           sendtmptype->is_replayable() ? sendtmpcount : sendtmpcount * sendtmptype->size(),
+                           recvtype->is_replayable() ? recvcount : recvcount * recvtype->size(),
+                           simgrid::smpi::Datatype::encode(sendtmptype), simgrid::smpi::Datatype::encode(recvtype)));
 
     retval = simgrid::smpi::Colls::alltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm);
 
@@ -610,7 +619,7 @@ int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype
              recvdisps == nullptr) {
     retval = MPI_ERR_ARG;
   } else {
-    int rank                           = simgrid::s4u::Actor::self()->getPid();
+    int rank                           = simgrid::s4u::this_actor::get_pid();
     int size               = comm->size();
     int send_size                      = 0;
     int recv_size                      = 0;
@@ -647,9 +656,10 @@ int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype
       trace_sendcounts->push_back(sendtmpcounts[i] * dt_size_send);
     }
 
-    TRACE_smpi_comm_in(rank, __FUNCTION__, new simgrid::instr::VarCollTIData(
-                                               "allToAllV", -1, send_size, trace_sendcounts, recv_size,
-                                               trace_recvcounts, encode_datatype(sendtype), encode_datatype(recvtype)));
+    TRACE_smpi_comm_in(rank, __func__,
+                       new simgrid::instr::VarCollTIData("allToAllV", -1, send_size, trace_sendcounts, recv_size,
+                                                         trace_recvcounts, simgrid::smpi::Datatype::encode(sendtype),
+                                                         simgrid::smpi::Datatype::encode(recvtype)));
 
     retval = simgrid::smpi::Colls::alltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
                                     recvdisps, recvtype, comm);
@@ -665,5 +675,3 @@ int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype
   smpi_bench_begin();
   return retval;
 }
-
-}