Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Useless test; TODO--.
[simgrid.git] / src / smpi / bindings / smpi_f77_coll.cpp
index 761ed29..ef64d76 100644 (file)
@@ -1,4 +1,4 @@
-/* Copyright (c) 2010-2019. The SimGrid Team. All rights reserved.          */
+/* Copyright (c) 2010-2021. 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. */
@@ -105,7 +105,6 @@ void mpi_alltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendtyp
 }
 
 void mpi_reduce_local_ (void *inbuf, void *inoutbuf, int* count, int* datatype, int* op, int* ierr){
-
  *ierr = MPI_Reduce_local(inbuf, inoutbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op));
 }
 
@@ -120,18 +119,16 @@ void mpi_reduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, in
 void mpi_alltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* old_sendtypes, void *recvbuf, int *recvcnts,
                       int *rdispls, int* old_recvtypes, int* comm, int* ierr){
   int size = simgrid::smpi::Comm::f2c(*comm)->size();
-  MPI_Datatype* sendtypes = new MPI_Datatype[size];
-  MPI_Datatype* recvtypes = new MPI_Datatype[size];
+  std::vector<MPI_Datatype> sendtypes(size);
+  std::vector<MPI_Datatype> recvtypes(size);
   for(int i=0; i< size; i++){
     if(FORT_IN_PLACE(sendbuf)!=MPI_IN_PLACE)
       sendtypes[i] = simgrid::smpi::Datatype::f2c(old_sendtypes[i]);
     recvtypes[i] = simgrid::smpi::Datatype::f2c(old_recvtypes[i]);
   }
   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
- *ierr = MPI_Alltoallw( sendbuf, sendcnts, sdispls, sendtypes, recvbuf, recvcnts, rdispls,
-                        recvtypes, simgrid::smpi::Comm::f2c(*comm));
-  delete[] sendtypes;
-  delete[] recvtypes;
+  *ierr   = MPI_Alltoallw(sendbuf, sendcnts, sdispls, sendtypes.data(), recvbuf, recvcnts, rdispls, recvtypes.data(),
+                        simgrid::smpi::Comm::f2c(*comm));
 }
 
 void mpi_exscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr){
@@ -142,7 +139,7 @@ void mpi_ibarrier_(int* comm, int* request, int* ierr) {
   MPI_Request req;
   *ierr = MPI_Ibarrier(simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -150,7 +147,7 @@ void mpi_ibcast_(void *buf, int* count, int* datatype, int* root, int* comm, int
   MPI_Request req;
   *ierr = MPI_Ibcast(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -161,7 +158,7 @@ void mpi_ireduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int*
   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
   *ierr = MPI_Ireduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -170,7 +167,7 @@ void mpi_iallreduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, in
   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
   *ierr = MPI_Iallreduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -180,7 +177,7 @@ void mpi_ireduce_scatter_(void* sendbuf, void* recvbuf, int* recvcounts, int* da
   *ierr = MPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, simgrid::smpi::Datatype::f2c(*datatype),
                         simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -191,7 +188,7 @@ void mpi_iscatter_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf,
   *ierr = MPI_Iscatter(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -202,7 +199,7 @@ void mpi_iscatterv_(void* sendbuf, int* sendcounts, int* displs, int* sendtype,
   *ierr = MPI_Iscatterv(sendbuf, sendcounts, displs, simgrid::smpi::Datatype::f2c(*sendtype),
                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -215,7 +212,7 @@ void mpi_igather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, i
   *ierr = MPI_Igather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                      recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -228,7 +225,7 @@ void mpi_igatherv_(void* sendbuf, int* sendcount, int* sendtype,
   *ierr = MPI_Igatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                      recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -239,7 +236,7 @@ void mpi_iallgather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf
   *ierr = MPI_Iallgather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                         recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -250,7 +247,7 @@ void mpi_iallgatherv_(void* sendbuf, int* sendcount, int* sendtype,
   *ierr = MPI_Iallgatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                         recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -260,7 +257,7 @@ void mpi_iscan_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op
   *ierr = MPI_Iscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype),
                    simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -271,7 +268,7 @@ void mpi_ialltoall_(void* sendbuf, int* sendcount, int* sendtype,
   *ierr = MPI_Ialltoall(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
                        recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -282,7 +279,7 @@ void mpi_ialltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendty
   *ierr = MPI_Ialltoallv(sendbuf, sendcounts, senddisps, simgrid::smpi::Datatype::f2c(*sendtype),
                        recvbuf, recvcounts, recvdisps, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -294,7 +291,7 @@ void mpi_ireduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, i
  *ierr = MPI_Ireduce_scatter_block(sendbuf, recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op),
                                   simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }
 
@@ -302,21 +299,19 @@ void mpi_ialltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* old_send
                       int *rdispls, int* old_recvtypes, int* comm, int* request, int* ierr){
   MPI_Request req;
   int size = simgrid::smpi::Comm::f2c(*comm)->size();
-  MPI_Datatype* sendtypes = new MPI_Datatype[size];
-  MPI_Datatype* recvtypes = new MPI_Datatype[size];
+  std::vector<MPI_Datatype> sendtypes(size);
+  std::vector<MPI_Datatype> recvtypes(size);
   for(int i=0; i< size; i++){
     if(FORT_IN_PLACE(sendbuf)!=MPI_IN_PLACE)
       sendtypes[i] = simgrid::smpi::Datatype::f2c(old_sendtypes[i]);
     recvtypes[i] = simgrid::smpi::Datatype::f2c(old_recvtypes[i]);
   }
   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
*ierr = MPI_Ialltoallw( sendbuf, sendcnts, sdispls, sendtypes, recvbuf, recvcnts, rdispls,
-                        recvtypes, simgrid::smpi::Comm::f2c(*comm), &req);
 *ierr   = MPI_Ialltoallw(sendbuf, sendcnts, sdispls, sendtypes.data(), recvbuf, recvcnts, rdispls, recvtypes.data(),
+                         simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
-  delete[] sendtypes;
-  delete[] recvtypes;
 }
 
 void mpi_iexscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr){
@@ -324,7 +319,7 @@ void mpi_iexscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int*
   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
  *ierr = MPI_Iexscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
   if(*ierr == MPI_SUCCESS) {
-    *request = req->add_f();
+    *request = req->c2f();
   }
 }