Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Have replay always use shared buffers instead of allocating new ones, even inside...
[simgrid.git] / src / smpi / colls / allreduce-rab1.c
index 0dc64c1..0a04e54 100644 (file)
@@ -1,4 +1,10 @@
-#include "colls.h"
+/* Copyright (c) 2013-2014. 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. */
+
+#include "colls_private.h"
 //#include <star-reduction.c>
 
 // NP pow of 2 for now
@@ -8,22 +14,24 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
 {
   MPI_Status status;
   MPI_Aint extent;
-  int tag = 4321, rank, nprocs, send_size, newcnt, share;
+  int tag = COLL_TAG_ALLREDUCE, rank, nprocs, send_size, newcnt, share;
   int pof2 = 1, mask, send_idx, recv_idx, dst, send_cnt, recv_cnt;
 
   void *recv, *tmp_buf;
 
-  MPI_Comm_rank(comm, &rank);
-  MPI_Comm_size(comm, &nprocs);
+  rank = smpi_comm_rank(comm);
+  nprocs = smpi_comm_size(comm);
+
+  if((nprocs&(nprocs-1)))
+    THROWF(arg_error,0, "allreduce rab1 algorithm can't be used with non power of two number of processes ! ");
 
-  MPI_Type_extent(dtype, &extent);
+  extent = smpi_datatype_get_extent(dtype);
 
   pof2 = 1;
   while (pof2 <= nprocs)
     pof2 <<= 1;
   pof2 >>= 1;
 
-  mask = 1;
   send_idx = recv_idx = 0;
 
   // uneven count
@@ -31,8 +39,8 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
     send_size = (count + nprocs) / nprocs;
     newcnt = send_size * nprocs;
 
-    recv = (void *) malloc(extent * newcnt);
-    tmp_buf = (void *) malloc(extent * newcnt);
+    recv = (void *) smpi_get_tmp_recvbuffer(extent * newcnt);
+    tmp_buf = (void *) smpi_get_tmp_sendbuffer(extent * newcnt);
     memcpy(recv, sbuff, extent * count);
 
 
@@ -47,10 +55,10 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
       else
         recv_idx = send_idx + (mask * share);
 
-      MPI_Sendrecv((char *) recv + send_idx * extent, send_cnt, dtype, dst, tag,
+      smpi_mpi_sendrecv((char *) recv + send_idx * extent, send_cnt, dtype, dst, tag,
                    tmp_buf, recv_cnt, dtype, dst, tag, comm, &status);
 
-      star_reduction(op, tmp_buf, (char *) recv + recv_idx * extent, &recv_cnt,
+      smpi_op_apply(op, tmp_buf, (char *) recv + recv_idx * extent, &recv_cnt,
                      &dtype);
 
       // update send_idx for next iteration 
@@ -59,16 +67,16 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
     }
 
     memcpy(tmp_buf, (char *) recv + recv_idx * extent, recv_cnt * extent);
-    MPI_Allgather(tmp_buf, recv_cnt, dtype, recv, recv_cnt, dtype, comm);
+    mpi_coll_allgather_fun(tmp_buf, recv_cnt, dtype, recv, recv_cnt, dtype, comm);
 
     memcpy(rbuff, recv, count * extent);
-    free(recv);
-    free(tmp_buf);
+    smpi_free_tmp_buffer(recv);
+    smpi_free_tmp_buffer(tmp_buf);
 
   }
 
   else {
-    tmp_buf = (void *) malloc(extent * count);
+    tmp_buf = (void *) smpi_get_tmp_sendbuffer(extent * count);
     memcpy(rbuff, sbuff, count * extent);
     mask = pof2 / 2;
     share = count / pof2;
@@ -81,10 +89,10 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
       else
         recv_idx = send_idx + (mask * share);
 
-      MPI_Sendrecv((char *) rbuff + send_idx * extent, send_cnt, dtype, dst,
+      smpi_mpi_sendrecv((char *) rbuff + send_idx * extent, send_cnt, dtype, dst,
                    tag, tmp_buf, recv_cnt, dtype, dst, tag, comm, &status);
 
-      star_reduction(op, tmp_buf, (char *) rbuff + recv_idx * extent, &recv_cnt,
+      smpi_op_apply(op, tmp_buf, (char *) rbuff + recv_idx * extent, &recv_cnt,
                      &dtype);
 
       // update send_idx for next iteration 
@@ -93,9 +101,9 @@ int smpi_coll_tuned_allreduce_rab1(void *sbuff, void *rbuff,
     }
 
     memcpy(tmp_buf, (char *) rbuff + recv_idx * extent, recv_cnt * extent);
-    MPI_Allgather(tmp_buf, recv_cnt, dtype, rbuff, recv_cnt, dtype, comm);
-    free(tmp_buf);
+    mpi_coll_allgather_fun(tmp_buf, recv_cnt, dtype, rbuff, recv_cnt, dtype, comm);
+    smpi_free_tmp_buffer(tmp_buf);
   }
 
-  return 0;
+  return MPI_SUCCESS;
 }