Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
switch way old SMP aware algos work, to be closer to the ones from mvapich
[simgrid.git] / src / smpi / colls / allreduce-smp-binomial.c
index 17dbd40..25f9837 100644 (file)
@@ -1,15 +1,16 @@
-#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"
 /* IMPLEMENTED BY PITCH PATARASUK 
    Non-topoloty-specific (however, number of cores/node need to be changed) 
    all-reduce operation designed for smp clusters
    It uses 2-layer communication: binomial for both intra-communication 
    inter-communication*/
 
-/* change number of core per smp-node
-   we assume that number of core per process will be the same for all implementations */
-#ifndef NUM_CORE
-#define NUM_CORE 8
-#endif
 
 /* ** NOTE **
    Use -DMPICH2 if this code does not compile.
@@ -32,26 +33,23 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
 {
   int comm_size, rank;
   void *tmp_buf;
-  int tag = 50;
+  int tag = COLL_TAG_ALLREDUCE;
   int mask, src, dst;
-  int num_core = NUM_CORE;
+
+  if(smpi_comm_get_leaders_comm(comm)==MPI_COMM_NULL){
+    smpi_comm_init_smp(comm);
+  }
+  int num_core=1;
+  if (smpi_comm_is_uniform(comm)){
+    num_core = smpi_comm_size(smpi_comm_get_intra_comm(comm));
+  }
   MPI_Status status;
-  /*
-     #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, &comm_size);
-  MPI_Comm_rank(comm, &rank);
-  MPI_Aint extent;
-  MPI_Type_extent(dtype, &extent);
-  tmp_buf = (void *) malloc(count * extent);
+
+  comm_size=smpi_comm_size(comm);
+  rank=smpi_comm_rank(comm);
+  MPI_Aint extent, lb;
+  smpi_datatype_extent(dtype, &lb, &extent);
+  tmp_buf = (void *) xbt_malloc(count * extent);
 
   /* compute intra and inter ranking */
   int intra_rank, inter_rank;
@@ -63,7 +61,7 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
   int inter_comm_size = (comm_size + num_core - 1) / num_core;
 
   /* copy input buffer to output buffer */
-  MPI_Sendrecv(send_buf, count, dtype, rank, tag,
+  smpi_mpi_sendrecv(send_buf, count, dtype, rank, tag,
                recv_buf, count, dtype, rank, tag, comm, &status);
 
   /* start binomial reduce intra communication inside each SMP node */
@@ -72,12 +70,12 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
     if ((mask & intra_rank) == 0) {
       src = (inter_rank * num_core) + (intra_rank | mask);
       if (src < comm_size) {
-        MPI_Recv(tmp_buf, count, dtype, src, tag, comm, &status);
-        star_reduction(op, tmp_buf, recv_buf, &count, &dtype);
+        smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
+        smpi_op_apply(op, tmp_buf, recv_buf, &count, &dtype);
       }
     } else {
       dst = (inter_rank * num_core) + (intra_rank & (~mask));
-      MPI_Send(recv_buf, count, dtype, dst, tag, comm);
+      smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
       break;
     }
     mask <<= 1;
@@ -91,12 +89,12 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
       if ((mask & inter_rank) == 0) {
         src = (inter_rank | mask) * num_core;
         if (src < comm_size) {
-          MPI_Recv(tmp_buf, count, dtype, src, tag, comm, &status);
-          star_reduction(op, tmp_buf, recv_buf, &count, &dtype);
+          smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
+          smpi_op_apply(op, tmp_buf, recv_buf, &count, &dtype);
         }
       } else {
         dst = (inter_rank & (~mask)) * num_core;
-        MPI_Send(recv_buf, count, dtype, dst, tag, comm);
+        smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
         break;
       }
       mask <<= 1;
@@ -110,7 +108,7 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
     while (mask < inter_comm_size) {
       if (inter_rank & mask) {
         src = (inter_rank - mask) * num_core;
-        MPI_Recv(recv_buf, count, dtype, src, tag, comm, &status);
+        smpi_mpi_recv(recv_buf, count, dtype, src, tag, comm, &status);
         break;
       }
       mask <<= 1;
@@ -121,7 +119,7 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
       if (inter_rank < inter_comm_size) {
         dst = (inter_rank + mask) * num_core;
         if (dst < comm_size) {
-          MPI_Send(recv_buf, count, dtype, dst, tag, comm);
+          smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
         }
       }
       mask >>= 1;
@@ -137,7 +135,7 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
   while (mask < num_core_in_current_smp) {
     if (intra_rank & mask) {
       src = (inter_rank * num_core) + (intra_rank - mask);
-      MPI_Recv(recv_buf, count, dtype, src, tag, comm, &status);
+      smpi_mpi_recv(recv_buf, count, dtype, src, tag, comm, &status);
       break;
     }
     mask <<= 1;
@@ -147,7 +145,7 @@ int smpi_coll_tuned_allreduce_smp_binomial(void *send_buf, void *recv_buf,
   while (mask > 0) {
     dst = (inter_rank * num_core) + (intra_rank + mask);
     if (dst < comm_size) {
-      MPI_Send(recv_buf, count, dtype, dst, tag, comm);
+      smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
     }
     mask >>= 1;
   }