Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add mpi alltoallv ring ring-light-barrier ring-mpi-barrier ring-one-barrier
authorPaul Bédaride <paul.bedaride@gmail.com>
Thu, 11 Apr 2013 14:20:05 +0000 (16:20 +0200)
committerPaul Bédaride <paul.bedaride@gmail.com>
Thu, 11 Apr 2013 14:20:05 +0000 (16:20 +0200)
buildtools/Cmake/AddTests.cmake
buildtools/Cmake/DefinePackages.cmake
src/smpi/colls/alltoallv-ring-light-barrier.c [new file with mode: 0644]
src/smpi/colls/alltoallv-ring-mpi-barrier.c [new file with mode: 0644]
src/smpi/colls/alltoallv-ring-one-barrier.c [new file with mode: 0644]
src/smpi/colls/alltoallv-ring.c [new file with mode: 0644]
src/smpi/colls/colls.h

index 0f99bb8..f72dd48 100644 (file)
@@ -384,7 +384,9 @@ if(NOT enable_memcheck)
                           simple bruck basic_linear pairwise)
         ADD_TEST(smpi-alltoall-coll-${ALLTOALL_COLL} ${CMAKE_BINARY_DIR}/bin/tesh ${TESH_OPTION} --cfg smpi/alltoall:${ALLTOALL_COLL} --cd ${CMAKE_BINARY_DIR}/teshsuite/smpi ${CMAKE_HOME_DIRECTORY}/teshsuite/smpi/alltoall_coll.tesh)
     ENDFOREACH()
-    FOREACH (ALLTOALLV_COLL default pair pair_light_barrier pair_mpi_barrier pair_one_barrier bruck pairwise)
+    FOREACH (ALLTOALLV_COLL default pair pair_light_barrier pair_mpi_barrier
+                           pair_one_barrier  ring ring_light_barrier
+                           ring_mpi_barrier ring_one_barrier bruck pairwise)
         ADD_TEST(smpi-alltoallv-coll-${ALLTOALLV_COLL} ${CMAKE_BINARY_DIR}/bin/tesh ${TESH_OPTION} --cfg smpi/alltoallv:${ALLTOALLV_COLL} --cd ${CMAKE_BINARY_DIR}/teshsuite/smpi ${CMAKE_HOME_DIRECTORY}/teshsuite/smpi/alltoallv_coll.tesh)
     ENDFOREACH()
     FOREACH (BCAST_COLL default arrival_nb arrival_pattern_aware arrival_pattern_aware_wait arrival_scatter
index f94813c..ff771aa 100644 (file)
@@ -158,6 +158,10 @@ set(SMPI_SRC
   src/smpi/colls/alltoallv-pair-light-barrier.c
   src/smpi/colls/alltoallv-pair-mpi-barrier.c
   src/smpi/colls/alltoallv-pair-one-barrier.c 
+  src/smpi/colls/alltoallv-ring.c
+  src/smpi/colls/alltoallv-ring-light-barrier.c
+  src/smpi/colls/alltoallv-ring-mpi-barrier.c
+  src/smpi/colls/alltoallv-ring-one-barrier.c
   src/smpi/colls/alltoallv-pairwise.c
   src/smpi/colls/alltoallv-bruck.c
   src/smpi/colls/bcast-arrival-nb.c
diff --git a/src/smpi/colls/alltoallv-ring-light-barrier.c b/src/smpi/colls/alltoallv-ring-light-barrier.c
new file mode 100644 (file)
index 0000000..7a42ab9
--- /dev/null
@@ -0,0 +1,67 @@
+#include "colls_private.h"
+/*****************************************************************************
+
+ * Function: alltoall_ring_light_barrier
+
+ * Return: int
+
+ * Inputs:
+    send_buff: send input buffer
+    send_count: number of elements to send
+    send_type: data type of elements being sent
+    recv_buff: receive output buffer
+    recv_count: number of elements to received
+    recv_type: data type of elements being received
+    comm: communicator
+
+ * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
+           Light barriers are inserted between communications in different
+           phases.
+
+ * Auther: Ahmad Faraj
+
+ ****************************************************************************/
+int
+smpi_coll_tuned_alltoallv_ring_light_barrier(void *send_buff, int *send_counts, int *send_disps,
+                                            MPI_Datatype send_type,
+                                            void *recv_buff, int *recv_counts, int *recv_disps,
+                                            MPI_Datatype recv_type,
+                                            MPI_Comm comm)
+{
+  MPI_Aint send_chunk, recv_chunk;
+  MPI_Status s;
+  int i, src, dst, rank, num_procs, next_dst, next_src;
+  int tag = 1;
+
+  char send_sync = 'a', recv_sync = 'b';
+  char *send_ptr = (char *) send_buff;
+  char *recv_ptr = (char *) recv_buff;
+
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  send_chunk = smpi_datatype_get_extent(send_type);
+  recv_chunk = smpi_datatype_get_extent(recv_type);
+
+  smpi_mpi_sendrecv(send_ptr + send_disps[rank] * send_chunk, send_counts[rank], send_type, rank, tag,
+               recv_ptr + recv_disps[rank] * recv_chunk, recv_counts[rank], recv_type, rank, tag,
+               comm, &s);
+
+  for (i = 1; i < num_procs; i++) {
+    src = (rank - i + num_procs) % num_procs;
+    dst = (rank + i) % num_procs;
+
+    smpi_mpi_sendrecv(send_ptr + send_disps[dst] * send_chunk, send_counts[dst], send_type,
+                 dst, tag, recv_ptr + recv_disps[src] * recv_chunk, recv_counts[src],
+                 recv_type, src, tag, comm, &s);
+
+    if ((i + 1) < num_procs) {
+      next_src = (rank - (i + 1) + num_procs) % num_procs;
+      next_dst = (rank + (i + 1) + num_procs) % num_procs;
+      smpi_mpi_sendrecv(&send_sync, 1, MPI_CHAR, next_src, tag,
+                   &recv_sync, 1, MPI_CHAR, next_dst, tag, comm, &s);
+
+    }
+  }
+
+  return MPI_SUCCESS;
+}
diff --git a/src/smpi/colls/alltoallv-ring-mpi-barrier.c b/src/smpi/colls/alltoallv-ring-mpi-barrier.c
new file mode 100644 (file)
index 0000000..273f15c
--- /dev/null
@@ -0,0 +1,53 @@
+#include "colls_private.h"
+/*****************************************************************************
+
+ * Function: alltoall_ring_mpi_barrier
+
+ * Return: int
+
+ * Inputs:
+    send_buff: send input buffer
+    send_count: number of elements to send
+    send_type: data type of elements being sent
+    recv_buff: receive output buffer
+    recv_count: number of elements to received
+    recv_type: data type of elements being received
+    comm: communicator
+
+ * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
+           MPI barriers are added between each two phases.
+
+ * Auther: Ahmad Faraj
+
+ ****************************************************************************/
+int
+smpi_coll_tuned_alltoallv_ring_mpi_barrier(void *send_buff, int *send_counts, int *send_disps,
+                                          MPI_Datatype send_type,
+                                          void *recv_buff, int *recv_counts, int *recv_disps,
+                                          MPI_Datatype recv_type, MPI_Comm comm)
+{
+  MPI_Status s;
+  MPI_Aint send_chunk, recv_chunk;
+  int i, src, dst, rank, num_procs;
+  int tag = 1;
+
+  char *send_ptr = (char *) send_buff;
+  char *recv_ptr = (char *) recv_buff;
+
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  send_chunk = smpi_datatype_get_extent(send_type);
+  recv_chunk = smpi_datatype_get_extent(recv_type);
+
+  for (i = 0; i < num_procs; i++) {
+    src = (rank - i + num_procs) % num_procs;
+    dst = (rank + i) % num_procs;
+
+    smpi_mpi_barrier(comm);
+    smpi_mpi_sendrecv(send_ptr + send_disps[dst] * send_chunk, send_counts[dst], send_type, dst,
+                 tag, recv_ptr + recv_disps[src] * recv_chunk, recv_counts[src], recv_type,
+                 src, tag, comm, &s);
+  }
+
+  return MPI_SUCCESS;
+}
diff --git a/src/smpi/colls/alltoallv-ring-one-barrier.c b/src/smpi/colls/alltoallv-ring-one-barrier.c
new file mode 100644 (file)
index 0000000..8e2d0cd
--- /dev/null
@@ -0,0 +1,51 @@
+#include "colls_private.h"
+/*****************************************************************************
+
+ * Function: alltoall_ring
+
+ * Return: int
+
+ * Inputs:
+    send_buff: send input buffer
+    send_count: number of elements to send
+    send_type: data type of elements being sent
+    recv_buff: receive output buffer
+    recv_count: number of elements to received
+    recv_type: data type of elements being received
+    comm: communicator
+
+ * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
+
+ * Auther: Ahmad Faraj
+
+ ****************************************************************************/
+int
+smpi_coll_tuned_alltoallv_ring_one_barrier(void *send_buff, int *send_counts, int *send_disps,
+                                          MPI_Datatype send_type,
+                                          void *recv_buff, int *recv_counts, int *recv_disps,
+                                          MPI_Datatype recv_type, MPI_Comm comm)
+{
+  MPI_Status s;
+  MPI_Aint send_chunk, recv_chunk;
+  int i, src, dst, rank, num_procs;
+  int tag = 1;
+
+  char *send_ptr = (char *) send_buff;
+  char *recv_ptr = (char *) recv_buff;
+
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  send_chunk = smpi_datatype_get_extent(send_type);
+  recv_chunk = smpi_datatype_get_extent(recv_type);
+
+  smpi_mpi_barrier(comm);
+  for (i = 0; i < num_procs; i++) {
+    src = (rank - i + num_procs) % num_procs;
+    dst = (rank + i) % num_procs;
+
+    smpi_mpi_sendrecv(send_ptr + send_disps[dst] * send_chunk, send_counts[dst], send_type, dst,
+                 tag, recv_ptr + recv_disps[src] * recv_chunk, recv_counts[src], recv_type,
+                 src, tag, comm, &s);
+  }
+  return MPI_SUCCESS;
+}
diff --git a/src/smpi/colls/alltoallv-ring.c b/src/smpi/colls/alltoallv-ring.c
new file mode 100644 (file)
index 0000000..e958bb4
--- /dev/null
@@ -0,0 +1,51 @@
+#include "colls_private.h"
+/*****************************************************************************
+
+ * Function: alltoall_ring
+
+ * Return: int
+
+ * Inputs:
+    send_buff: send input buffer
+    send_count: number of elements to send
+    send_type: data type of elements being sent
+    recv_buff: receive output buffer
+    recv_count: number of elements to received
+    recv_type: data type of elements being received
+    comm: communicator
+
+ * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
+
+ * Auther: Ahmad Faraj
+
+ ****************************************************************************/
+int
+smpi_coll_tuned_alltoallv_ring(void *send_buff, int *send_counts, int *send_disps,
+                             MPI_Datatype send_type,
+                             void *recv_buff,int *recv_counts, int *recv_disps, 
+                             MPI_Datatype recv_type,
+                              MPI_Comm comm)
+{
+  MPI_Status s;
+  MPI_Aint send_chunk, recv_chunk;
+  int i, src, dst, rank, num_procs;
+  int tag = 11;
+
+  char *send_ptr = (char *) send_buff;
+  char *recv_ptr = (char *) recv_buff;
+
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  send_chunk = smpi_datatype_get_extent(send_type);
+  recv_chunk = smpi_datatype_get_extent(recv_type);
+
+  for (i = 0; i < num_procs; i++) {
+    src = (rank - i + num_procs) % num_procs;
+    dst = (rank + i) % num_procs;
+
+    smpi_mpi_sendrecv(send_ptr + send_disps[dst] * send_chunk, send_counts[dst], send_type, dst,
+                 tag, recv_ptr + recv_disps[src] * recv_chunk, recv_counts[src], recv_type,
+                 src, tag, comm, &s);
+  }
+  return MPI_SUCCESS;
+}
index 043f15f..8349bc5 100644 (file)
@@ -115,6 +115,10 @@ COLL_APPLY(action, COLL_ALLTOALLV_SIG, pair) COLL_sep \
 COLL_APPLY(action, COLL_ALLTOALLV_SIG, pair_light_barrier) COLL_sep \
 COLL_APPLY(action, COLL_ALLTOALLV_SIG, pair_mpi_barrier) COLL_sep \
 COLL_APPLY(action, COLL_ALLTOALLV_SIG, pair_one_barrier) COLL_sep \
+COLL_APPLY(action, COLL_ALLTOALLV_SIG, ring) COLL_sep \
+COLL_APPLY(action, COLL_ALLTOALLV_SIG, ring_light_barrier) COLL_sep \
+COLL_APPLY(action, COLL_ALLTOALLV_SIG, ring_mpi_barrier) COLL_sep \
+COLL_APPLY(action, COLL_ALLTOALLV_SIG, ring_one_barrier) COLL_sep \
 COLL_APPLY(action, COLL_ALLTOALLV_SIG, pairwise)
 
 COLL_ALLTOALLVS(COLL_PROTO, COLL_NOsep)