Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add scatter algos from ompi
[simgrid.git] / src / smpi / colls / allreduce-rdb.c
index 5e3cf46..44717ad 100644 (file)
@@ -1,12 +1,12 @@
-#include "colls.h"
+#include "colls_private.h"
 //#include <star-reduction.c>
 
 int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
                                   MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)
 {
-  int nprocs, rank, type_size, tag = 543;
+  int nprocs, rank, tag = 543;
   int mask, dst, pof2, newrank, rem, newdst;
-  MPI_Aint extent;
+  MPI_Aint extent, lb;
   MPI_Status status;
   void *tmp_buf = NULL;
   /*
@@ -19,20 +19,14 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
      uop  = op_ptr->op;
      #endif
    */
-  MPI_Comm_size(comm, &nprocs);
-  MPI_Comm_rank(comm, &rank);
-
-  MPI_Type_extent(dtype, &extent);
-  tmp_buf = (void *) malloc(count * extent);
-  if (!tmp_buf) {
-    printf("Could not allocate memory for tmp_buf\n");
-    return 1;
-  }
+  nprocs=smpi_comm_size(comm);
+  rank=smpi_comm_rank(comm);
 
-  MPI_Sendrecv(sbuff, count, dtype, rank, 500,
-               rbuff, count, dtype, rank, 500, comm, &status);
+  smpi_datatype_extent(dtype, &lb, &extent);
+  tmp_buf = (void *) xbt_malloc(count * extent);
 
-  MPI_Type_size(dtype, &type_size);
+  smpi_mpi_sendrecv(sbuff, count, dtype, rank, 500,
+               rbuff, count, dtype, rank, 500, comm, &status);
 
   // find nearest power-of-two less than or equal to comm_size
   pof2 = 1;
@@ -52,7 +46,7 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
     // even       
     if (rank % 2 == 0) {
 
-      MPI_Send(rbuff, count, dtype, rank + 1, tag, comm);
+      smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
 
       // temporarily set the rank to -1 so that this
       // process does not pariticipate in recursive
@@ -60,11 +54,11 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
       newrank = -1;
     } else                      // odd
     {
-      MPI_Recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
+      smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
       // do the reduction on received data. since the
       // ordering is right, it doesn't matter whether
       // the operation is commutative or not.
-      star_reduction(op, tmp_buf, rbuff, &count, &dtype);
+      smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
 
       // change the rank 
       newrank = rank / 2;
@@ -92,7 +86,7 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
 
       // Send the most current data, which is in recvbuf. Recv
       // into tmp_buf 
-      MPI_Sendrecv(rbuff, count, dtype, dst, tag, tmp_buf, count, dtype,
+      smpi_mpi_sendrecv(rbuff, count, dtype, dst, tag, tmp_buf, count, dtype,
                    dst, tag, comm, &status);
 
       // tmp_buf contains data received in this step.
@@ -102,13 +96,13 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
       // we assume it is commuttive op
       //      if (op -> op_commute  || (dst < rank))
       if ((dst < rank)) {
-        star_reduction(op, tmp_buf, rbuff, &count, &dtype);
+        smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
       } else                    // op is noncommutative and the order is not right
       {
-        star_reduction(op, rbuff, tmp_buf, &count, &dtype);
+        smpi_op_apply(op, rbuff, tmp_buf, &count, &dtype);
 
         // copy result back into recvbuf
-        MPI_Sendrecv(tmp_buf, count, dtype, rank, tag, rbuff, count,
+        smpi_mpi_sendrecv(tmp_buf, count, dtype, rank, tag, rbuff, count,
                      dtype, rank, tag, comm, &status);
       }
       mask <<= 1;
@@ -120,11 +114,11 @@ int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
 
   if (rank < 2 * rem) {
     if (rank % 2)               // odd 
-      MPI_Send(rbuff, count, dtype, rank - 1, tag, comm);
+      smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
     else                        // even 
-      MPI_Recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
+      smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
   }
 
   free(tmp_buf);
-  return 0;
+  return MPI_SUCCESS;
 }