Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Replace dynamic C-style arrays with std::vector.
authorArnaud Giersch <arnaud.giersch@univ-fcomte.fr>
Fri, 20 Nov 2020 20:54:06 +0000 (21:54 +0100)
committerArnaud Giersch <arnaud.giersch@univ-fcomte.fr>
Sat, 21 Nov 2020 17:01:37 +0000 (18:01 +0100)
src/smpi/bindings/smpi_f77_request.cpp
src/smpi/bindings/smpi_f77_type.cpp
src/smpi/bindings/smpi_pmpi_request.cpp
src/smpi/mpi/smpi_comm.cpp
src/smpi/mpi/smpi_win.cpp

index 0a1cf7d..cfead1e 100644 (file)
@@ -150,15 +150,13 @@ void mpi_start_(int* request, int* ierr) {
   *ierr = MPI_Start(&req);
 }
 
-void mpi_startall_(int* count, int* requests, int* ierr) {
-  MPI_Request* reqs;
-
-  reqs = xbt_new(MPI_Request, *count);
+void mpi_startall_(int* count, int* requests, int* ierr)
+{
+  std::vector<MPI_Request> reqs(*count);
   for (int i = 0; i < *count; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr = MPI_Startall(*count, reqs);
-  xbt_free(reqs);
+  *ierr = MPI_Startall(*count, reqs.data());
 }
 
 void mpi_wait_(int* request, MPI_Status* status, int* ierr) {
@@ -171,14 +169,13 @@ void mpi_wait_(int* request, MPI_Status* status, int* ierr) {
    }
 }
 
-void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int* ierr) {
-  MPI_Request* reqs;
-
-  reqs = xbt_new(MPI_Request, *count);
+void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int* ierr)
+{
+  std::vector<MPI_Request> reqs(*count);
   for (int i = 0; i < *count; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr = MPI_Waitany(*count, reqs, index, status);
+  *ierr = MPI_Waitany(*count, reqs.data(), index, status);
   if(*index!=MPI_UNDEFINED){
     if(reqs[*index]==MPI_REQUEST_NULL){
         simgrid::smpi::Request::free_f(requests[*index]);
@@ -186,46 +183,37 @@ void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int
     }
   *index=*index+1;
   }
-  xbt_free(reqs);
 }
 
-void mpi_waitall_(int* count, int* requests, MPI_Status* status, int* ierr) {
-  MPI_Request* reqs;
-  int i;
-
-  reqs = xbt_new(MPI_Request, *count);
-  for(i = 0; i < *count; i++) {
+void mpi_waitall_(int* count, int* requests, MPI_Status* status, int* ierr)
+{
+  std::vector<MPI_Request> reqs(*count);
+  for (int i = 0; i < *count; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr = MPI_Waitall(*count, reqs, FORT_STATUSES_IGNORE(status));
-  for(i = 0; i < *count; i++) {
-      if(reqs[i]==MPI_REQUEST_NULL){
-          simgrid::smpi::Request::free_f(requests[i]);
-          requests[i]=MPI_FORTRAN_REQUEST_NULL;
-      }
+  *ierr = MPI_Waitall(*count, reqs.data(), FORT_STATUSES_IGNORE(status));
+  for (int i = 0; i < *count; i++) {
+    if (reqs[i] == MPI_REQUEST_NULL) {
+      simgrid::smpi::Request::free_f(requests[i]);
+      requests[i] = MPI_FORTRAN_REQUEST_NULL;
+    }
   }
-
-  xbt_free(reqs);
 }
 
 void mpi_waitsome_ (int* incount, int* requests, int *outcount, int *indices, MPI_Status* status, int* ierr)
 {
-  MPI_Request* reqs;
-  int i;
-
-  reqs = xbt_new(MPI_Request, *incount);
-  for(i = 0; i < *incount; i++) {
+  std::vector<MPI_Request> reqs(*incount);
+  for (int i = 0; i < *incount; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr = MPI_Waitsome(*incount, reqs, outcount, indices, status);
-  for(i=0;i<*outcount;i++){
+  *ierr = MPI_Waitsome(*incount, reqs.data(), outcount, indices, status);
+  for (int i = 0; i < *outcount; i++) {
     if(reqs[indices[i]]==MPI_REQUEST_NULL){
         simgrid::smpi::Request::free_f(requests[indices[i]]);
         requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
     }
     indices[i]++;
   }
-  xbt_free(reqs);
 }
 
 void mpi_test_ (int * request, int *flag, MPI_Status * status, int* ierr){
@@ -237,31 +225,28 @@ void mpi_test_ (int * request, int *flag, MPI_Status * status, int* ierr){
   }
 }
 
-void mpi_testall_ (int* count, int * requests,  int *flag, MPI_Status * statuses, int* ierr){
-  int i;
-  MPI_Request* reqs = xbt_new(MPI_Request, *count);
-  for(i = 0; i < *count; i++) {
+void mpi_testall_(int* count, int* requests, int* flag, MPI_Status* statuses, int* ierr)
+{
+  std::vector<MPI_Request> reqs(*count);
+  for (int i = 0; i < *count; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr= MPI_Testall(*count, reqs, flag, FORT_STATUSES_IGNORE(statuses));
-  for(i = 0; i < *count; i++) {
+  *ierr = MPI_Testall(*count, reqs.data(), flag, FORT_STATUSES_IGNORE(statuses));
+  for (int i = 0; i < *count; i++) {
     if(reqs[i]==MPI_REQUEST_NULL){
         simgrid::smpi::Request::free_f(requests[i]);
         requests[i]=MPI_FORTRAN_REQUEST_NULL;
     }
   }
-  xbt_free(reqs);
 }
 
 void mpi_testany_ (int* count, int* requests, int *index, int *flag, MPI_Status* status, int* ierr)
 {
-  MPI_Request* reqs;
-
-  reqs = xbt_new(MPI_Request, *count);
+  std::vector<MPI_Request> reqs(*count);
   for (int i = 0; i < *count; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
   }
-  *ierr = MPI_Testany(*count, reqs, index, flag, FORT_STATUS_IGNORE(status));
+  *ierr = MPI_Testany(*count, reqs.data(), index, flag, FORT_STATUS_IGNORE(status));
   if(*index!=MPI_UNDEFINED){
     if(reqs[*index]==MPI_REQUEST_NULL){
     simgrid::smpi::Request::free_f(requests[*index]);
@@ -269,27 +254,23 @@ void mpi_testany_ (int* count, int* requests, int *index, int *flag, MPI_Status*
     }
   *index=*index+1;
   }
-  xbt_free(reqs);
 }
 
-void mpi_testsome_ (int* incount, int*  requests, int* outcount, int* indices, MPI_Status*  statuses, int* ierr) {
-  MPI_Request* reqs;
-  int i;
-
-  reqs = xbt_new(MPI_Request, *incount);
-  for(i = 0; i < *incount; i++) {
+void mpi_testsome_(int* incount, int* requests, int* outcount, int* indices, MPI_Status* statuses, int* ierr)
+{
+  std::vector<MPI_Request> reqs(*incount);
+  for (int i = 0; i < *incount; i++) {
     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
     indices[i]=0;
   }
-  *ierr = MPI_Testsome(*incount, reqs, outcount, indices, FORT_STATUSES_IGNORE(statuses));
-  for(i=0;i<*incount;i++){
+  *ierr = MPI_Testsome(*incount, reqs.data(), outcount, indices, FORT_STATUSES_IGNORE(statuses));
+  for (int i = 0; i < *incount; i++) {
     if(reqs[indices[i]]==MPI_REQUEST_NULL){
       simgrid::smpi::Request::free_f(requests[indices[i]]);
       requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
     }
     indices[i]++;
   }
-  xbt_free(reqs);
 }
 
 void mpi_probe_ (int* source, int* tag, int* comm, MPI_Status*  status, int* ierr) {
index b492a10..648c7d5 100644 (file)
@@ -186,31 +186,28 @@ void mpi_type_create_indexed_block_ (int* count, int* blocklength, int* indices,
 
 void mpi_type_struct_ (int* count, int* blocklens, int* indices, int* old_types, int*  newtype, int* ierr) {
   MPI_Datatype tmp;
-  auto* types        = static_cast<MPI_Datatype*>(xbt_malloc(*count * sizeof(MPI_Datatype)));
-  auto* indices_aint = new MPI_Aint[*count];
+  std::vector<MPI_Aint> indices_aint(*count);
+  std::vector<MPI_Datatype> types(*count);
   for (int i = 0; i < *count; i++) {
     indices_aint[i]=indices[i];
     types[i] = simgrid::smpi::Datatype::f2c(old_types[i]);
   }
-  *ierr = MPI_Type_struct(*count, blocklens, indices_aint, types, &tmp);
+  *ierr = MPI_Type_struct(*count, blocklens, indices_aint.data(), types.data(), &tmp);
   if(*ierr == MPI_SUCCESS) {
     *newtype = tmp->add_f();
   }
-  xbt_free(types);
-  delete[] indices_aint;
 }
 
 void mpi_type_create_struct_(int* count, int* blocklens, MPI_Aint* indices, int*  old_types, int*  newtype, int* ierr){
   MPI_Datatype tmp;
-  auto* types = static_cast<MPI_Datatype*>(xbt_malloc(*count * sizeof(MPI_Datatype)));
+  std::vector<MPI_Datatype> types(*count);
   for (int i = 0; i < *count; i++) {
     types[i] = simgrid::smpi::Datatype::f2c(old_types[i]);
   }
-  *ierr = MPI_Type_create_struct(*count, blocklens, indices, types, &tmp);
+  *ierr = MPI_Type_create_struct(*count, blocklens, indices, types.data(), &tmp);
   if(*ierr == MPI_SUCCESS) {
     *newtype = tmp->add_f();
   }
-  xbt_free(types);
 }
 
 void mpi_pack_ (void* inbuf, int* incount, int* type, void* outbuf, int* outcount, int* position, int* comm, int* ierr) {
index 7d6a45d..78db9c0 100644 (file)
@@ -461,12 +461,12 @@ int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst,
 
   int size = datatype->get_extent() * count;
   xbt_assert(size > 0);
-  void* recvbuf = xbt_new0(char, size);
-  retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
+  std::vector<char> recvbuf(size);
+  retval =
+      MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf.data(), count, datatype, src, recvtag, comm, status);
   if(retval==MPI_SUCCESS){
-    simgrid::smpi::Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
+    simgrid::smpi::Datatype::copy(recvbuf.data(), count, datatype, buf, count, datatype);
   }
-  xbt_free(recvbuf);
   return retval;
 }
 
index 307fcdd..e1cf7e6 100644 (file)
@@ -233,7 +233,6 @@ MPI_Comm Comm::split(int color, int key)
   if (this == MPI_COMM_UNINITIALIZED)
     return smpi_process()->comm_world()->split(color, key);
   int system_tag = -123;
-  int* recvbuf;
 
   MPI_Group group_root = nullptr;
   MPI_Group group_out  = nullptr;
@@ -242,15 +241,13 @@ MPI_Comm Comm::split(int color, int key)
   int size             = this->size();
   /* Gather all colors and keys on rank 0 */
   const std::array<int, 2> sendbuf = {{color, key}};
-  if (myrank == 0) {
-    recvbuf = xbt_new(int, 2 * size);
-  } else {
-    recvbuf = nullptr;
-  }
-  gather__default(sendbuf.data(), 2, MPI_INT, recvbuf, 2, MPI_INT, 0, this);
+  std::vector<int> recvbuf;
+  if (myrank == 0)
+    recvbuf.resize(2 * size);
+  gather__default(sendbuf.data(), 2, MPI_INT, recvbuf.data(), 2, MPI_INT, 0, this);
   /* Do the actual job */
   if (myrank == 0) {
-    MPI_Group* group_snd = xbt_new(MPI_Group, size);
+    std::vector<MPI_Group> group_snd(size);
     std::vector<std::pair<int, int>> rankmap;
     rankmap.reserve(size);
     for (int i = 0; i < size; i++) {
@@ -274,7 +271,7 @@ MPI_Comm Comm::split(int color, int key)
           s4u::Actor* actor = group->actor(rankmap[j].second);
           group_out->set_mapping(actor, j);
         }
-        MPI_Request* requests = xbt_new(MPI_Request, rankmap.size());
+        std::vector<MPI_Request> requests(rankmap.size());
         int reqs              = 0;
         for (auto const& rank : rankmap) {
           if (rank.second != 0) {
@@ -286,12 +283,9 @@ MPI_Comm Comm::split(int color, int key)
         if(i != 0 && group_out != MPI_COMM_WORLD->group() && group_out != MPI_GROUP_EMPTY)
           Group::unref(group_out);
 
-        Request::waitall(reqs, requests, MPI_STATUS_IGNORE);
-        xbt_free(requests);
+        Request::waitall(reqs, requests.data(), MPI_STATUS_IGNORE);
       }
     }
-    xbt_free(recvbuf);
-    xbt_free(group_snd);
     group_out = group_root; /* exit with root's group */
   } else {
     if(color != MPI_UNDEFINED) {
index 806b363..6f1245e 100644 (file)
@@ -469,7 +469,7 @@ int Win::start(MPI_Group group, int /*assert*/)
   int i             = 0;
   int j             = 0;
   int size          = group->size();
-  MPI_Request* reqs = xbt_new0(MPI_Request, size);
+  std::vector<MPI_Request> reqs(size);
 
   XBT_DEBUG("Entering MPI_Win_Start");
   while (j != size) {
@@ -481,12 +481,11 @@ int Win::start(MPI_Group group, int /*assert*/)
     j++;
   }
   size = i;
-  Request::startall(size, reqs);
-  Request::waitall(size, reqs, MPI_STATUSES_IGNORE);
+  Request::startall(size, reqs.data());
+  Request::waitall(size, reqs.data(), MPI_STATUSES_IGNORE);
   for (i = 0; i < size; i++) {
     Request::unref(&reqs[i]);
   }
-  xbt_free(reqs);
   opened_++; //we're open for business !
   group_=group;
   group->ref();
@@ -500,7 +499,7 @@ int Win::post(MPI_Group group, int /*assert*/)
   int i             = 0;
   int j             = 0;
   int size = group->size();
-  MPI_Request* reqs = xbt_new0(MPI_Request, size);
+  std::vector<MPI_Request> reqs(size);
 
   XBT_DEBUG("Entering MPI_Win_Post");
   while(j!=size){
@@ -513,12 +512,11 @@ int Win::post(MPI_Group group, int /*assert*/)
   }
   size=i;
 
-  Request::startall(size, reqs);
-  Request::waitall(size, reqs, MPI_STATUSES_IGNORE);
+  Request::startall(size, reqs.data());
+  Request::waitall(size, reqs.data(), MPI_STATUSES_IGNORE);
   for(i=0;i<size;i++){
     Request::unref(&reqs[i]);
   }
-  xbt_free(reqs);
   opened_++; //we're open for business !
   group_=group;
   group->ref();
@@ -534,7 +532,7 @@ int Win::complete(){
   int i             = 0;
   int j             = 0;
   int size          = group_->size();
-  MPI_Request* reqs = xbt_new0(MPI_Request, size);
+  std::vector<MPI_Request> reqs(size);
 
   while(j!=size){
     int dst = comm_->group()->rank(group_->actor(j));
@@ -546,13 +544,12 @@ int Win::complete(){
   }
   size=i;
   XBT_DEBUG("Win_complete - Sending sync messages to %d processes", size);
-  Request::startall(size, reqs);
-  Request::waitall(size, reqs, MPI_STATUSES_IGNORE);
+  Request::startall(size, reqs.data());
+  Request::waitall(size, reqs.data(), MPI_STATUSES_IGNORE);
 
   for(i=0;i<size;i++){
     Request::unref(&reqs[i]);
   }
-  xbt_free(reqs);
 
   int finished = finish_comms();
   XBT_DEBUG("Win_complete - Finished %d RMA calls", finished);
@@ -568,7 +565,7 @@ int Win::wait(){
   int i             = 0;
   int j             = 0;
   int size          = group_->size();
-  MPI_Request* reqs = xbt_new0(MPI_Request, size);
+  std::vector<MPI_Request> reqs(size);
 
   while(j!=size){
     int src = comm_->group()->rank(group_->actor(j));
@@ -580,12 +577,11 @@ int Win::wait(){
   }
   size=i;
   XBT_DEBUG("Win_wait - Receiving sync messages from %d processes", size);
-  Request::startall(size, reqs);
-  Request::waitall(size, reqs, MPI_STATUSES_IGNORE);
+  Request::startall(size, reqs.data());
+  Request::waitall(size, reqs.data(), MPI_STATUSES_IGNORE);
   for(i=0;i<size;i++){
     Request::unref(&reqs[i]);
   }
-  xbt_free(reqs);
   int finished = finish_comms();
   XBT_DEBUG("Win_wait - Finished %d RMA calls", finished);