Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add mpi alltoallv pair
authorPaul Bédaride <paul.bedaride@gmail.com>
Thu, 11 Apr 2013 12:56:21 +0000 (14:56 +0200)
committerPaul Bédaride <paul.bedaride@gmail.com>
Thu, 11 Apr 2013 12:56:21 +0000 (14:56 +0200)
src/smpi/colls/alltoallv-pair.c [new file with mode: 0644]

diff --git a/src/smpi/colls/alltoallv-pair.c b/src/smpi/colls/alltoallv-pair.c
new file mode 100644 (file)
index 0000000..50839df
--- /dev/null
@@ -0,0 +1,49 @@
+#include "colls_private.h"
+
+/*****************************************************************************
+
+ * Function: alltoall_pair
+
+ * 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 when P is power of two. In each phase of P - 1
+           phases, nodes in pair communicate their data.
+
+ * Auther: Ahmad Faraj
+
+ ****************************************************************************/
+int smpi_coll_tuned_alltoallv_pair(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;
+  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 = dst = rank ^ i;
+    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;
+}