Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add collectives for allgather, allreduce, bcast and reduce
[simgrid.git] / src / smpi / colls / allreduce-rdb.c
diff --git a/src/smpi/colls/allreduce-rdb.c b/src/smpi/colls/allreduce-rdb.c
new file mode 100644 (file)
index 0000000..5e3cf46
--- /dev/null
@@ -0,0 +1,130 @@
+#include "colls.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 mask, dst, pof2, newrank, rem, newdst;
+  MPI_Aint extent;
+  MPI_Status status;
+  void *tmp_buf = NULL;
+  /*
+     #ifdef MPICH2_REDUCTION
+     MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
+     #else
+     MPI_User_function *uop;
+     struct MPIR_OP *op_ptr;
+     op_ptr = MPIR_ToPointer(op);
+     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;
+  }
+
+  MPI_Sendrecv(sbuff, count, dtype, rank, 500,
+               rbuff, count, dtype, rank, 500, comm, &status);
+
+  MPI_Type_size(dtype, &type_size);
+
+  // find nearest power-of-two less than or equal to comm_size
+  pof2 = 1;
+  while (pof2 <= nprocs)
+    pof2 <<= 1;
+  pof2 >>= 1;
+
+  rem = nprocs - pof2;
+
+  // In the non-power-of-two case, all even-numbered
+  // processes of rank < 2*rem send their data to
+  // (rank+1). These even-numbered processes no longer
+  // participate in the algorithm until the very end. The
+  // remaining processes form a nice power-of-two. 
+
+  if (rank < 2 * rem) {
+    // even       
+    if (rank % 2 == 0) {
+
+      MPI_Send(rbuff, count, dtype, rank + 1, tag, comm);
+
+      // temporarily set the rank to -1 so that this
+      // process does not pariticipate in recursive
+      // doubling
+      newrank = -1;
+    } else                      // odd
+    {
+      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);
+
+      // change the rank 
+      newrank = rank / 2;
+    }
+  }
+
+  else                          // rank >= 2 * rem 
+    newrank = rank - rem;
+
+  // If op is user-defined or count is less than pof2, use
+  // recursive doubling algorithm. Otherwise do a reduce-scatter
+  // followed by allgather. (If op is user-defined,
+  // derived datatypes are allowed and the user could pass basic
+  // datatypes on one process and derived on another as long as
+  // the type maps are the same. Breaking up derived
+  // datatypes to do the reduce-scatter is tricky, therefore
+  // using recursive doubling in that case.) 
+
+  if (newrank != -1) {
+    mask = 0x1;
+    while (mask < pof2) {
+      newdst = newrank ^ mask;
+      // find real rank of dest 
+      dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
+
+      // Send the most current data, which is in recvbuf. Recv
+      // into tmp_buf 
+      MPI_Sendrecv(rbuff, count, dtype, dst, tag, tmp_buf, count, dtype,
+                   dst, tag, comm, &status);
+
+      // tmp_buf contains data received in this step.
+      // recvbuf contains data accumulated so far 
+
+      // op is commutative OR the order is already right
+      // we assume it is commuttive op
+      //      if (op -> op_commute  || (dst < rank))
+      if ((dst < rank)) {
+        star_reduction(op, tmp_buf, rbuff, &count, &dtype);
+      } else                    // op is noncommutative and the order is not right
+      {
+        star_reduction(op, rbuff, tmp_buf, &count, &dtype);
+
+        // copy result back into recvbuf
+        MPI_Sendrecv(tmp_buf, count, dtype, rank, tag, rbuff, count,
+                     dtype, rank, tag, comm, &status);
+      }
+      mask <<= 1;
+    }
+  }
+  // In the non-power-of-two case, all odd-numbered processes of
+  // rank < 2 * rem send the result to (rank-1), the ranks who didn't
+  // participate above.
+
+  if (rank < 2 * rem) {
+    if (rank % 2)               // odd 
+      MPI_Send(rbuff, count, dtype, rank - 1, tag, comm);
+    else                        // even 
+      MPI_Recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
+  }
+
+  free(tmp_buf);
+  return 0;
+}