Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add collectives for allgather, allreduce, bcast and reduce
[simgrid.git] / src / smpi / colls / reduce-arrival-pattern-aware.c
diff --git a/src/smpi/colls/reduce-arrival-pattern-aware.c b/src/smpi/colls/reduce-arrival-pattern-aware.c
new file mode 100644 (file)
index 0000000..0182a6b
--- /dev/null
@@ -0,0 +1,358 @@
+#include "colls.h"
+//#include <star-reduction.c>
+
+int reduce_arrival_pattern_aware_segment_size_in_byte = 8192;
+
+#ifndef HEADER_SIZE
+#define HEADER_SIZE 1024
+#endif
+
+#ifndef MAX_NODE
+#define MAX_NODE 1024
+#endif
+
+/* Non-topology-specific pipelined linear-reduce function */
+int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf,
+                                                 int count,
+                                                 MPI_Datatype datatype,
+                                                 MPI_Op op, int root,
+                                                 MPI_Comm comm)
+{
+  int rank;
+  MPI_Comm_rank(comm, &rank);
+
+  int tag = 50;
+  MPI_Status status;
+  MPI_Request request;
+  MPI_Request *send_request_array;
+  MPI_Request *recv_request_array;
+  MPI_Status *send_status_array;
+  MPI_Status *recv_status_array;
+
+  MPI_Status temp_status_array[MAX_NODE];
+
+  int size;
+  int i;
+
+  int sent_count;
+  int header_index;
+  int flag_array[MAX_NODE];
+  int already_received[MAX_NODE];
+
+  int header_buf[HEADER_SIZE];
+  char temp_buf[MAX_NODE];
+
+  MPI_Aint extent;
+  MPI_Type_extent(datatype, &extent);
+
+  /* source and destination */
+  int to, from;
+
+  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
+  MPI_Comm_size(MPI_COMM_WORLD, &size);
+
+
+  /* segment is segment size in number of elements (not bytes) */
+  int segment = reduce_arrival_pattern_aware_segment_size_in_byte / extent;
+
+  /* pipeline length */
+  int pipe_length = count / segment;
+
+  /* use for buffer offset for sending and receiving data = segment size in byte */
+  int increment = segment * extent;
+
+  /* if the input size is not divisible by segment size => 
+     the small remainder will be done with native implementation */
+  int remainder = count % segment;
+
+
+  /* value == 0 means root has not send data (or header) to the node yet */
+  for (i = 0; i < MAX_NODE; i++) {
+    already_received[i] = 0;
+  }
+
+  char *tmp_buf;
+  tmp_buf = (char *) malloc(count * extent);
+
+  MPI_Sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank,
+               tag, comm, &status);
+
+
+
+  /* when a message is smaller than a block size => no pipeline */
+  if (count <= segment) {
+
+    if (rank == 0) {
+      sent_count = 0;
+
+      while (sent_count < (size - 1)) {
+
+        for (i = 1; i < size; i++) {
+          if (already_received[i] == 0)
+            MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
+                       MPI_STATUSES_IGNORE);
+        }
+
+        header_index = 0;
+        /* recv 1-byte message */
+        for (i = 0; i < size; i++) {
+          if (i == rank)
+            continue;
+
+          /* 1-byte message arrive */
+          if ((flag_array[i] == 1) && (already_received[i] == 0)) {
+            MPI_Recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
+            header_buf[header_index] = i;
+            header_index++;
+            sent_count++;
+
+
+            //printf("root send to %d recv from %d : data = ",to,from);
+            /*
+               for (i=0;i<=header_index;i++) {
+               printf("%d ",header_buf[i]);
+               }
+               printf("\n");
+             */
+            /* will receive in the next step */
+            already_received[i] = 1;
+          }
+        }
+
+        /* send header followed by receive and reduce data */
+        if (header_index != 0) {
+          header_buf[header_index] = -1;
+          to = header_buf[0];
+          from = header_buf[header_index - 1];
+
+          MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
+          MPI_Recv(tmp_buf, count, datatype, from, tag, comm, &status);
+          star_reduction(op, tmp_buf, rbuf, &count, &datatype);
+        }
+      }                         /* while loop */
+    }
+
+    /* root */
+    /* non-root */
+    else {
+
+      /* send 1-byte message to root */
+      MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
+
+      /* wait for header and data, forward when required */
+      MPI_Recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
+               &status);
+      //      MPI_Recv(buf,count,datatype,MPI_ANY_SOURCE,tag,comm,&status);
+
+      /* search for where it is */
+      int myordering = 0;
+      while (rank != header_buf[myordering]) {
+        myordering++;
+      }
+
+      /* forward header */
+      if (header_buf[myordering + 1] != -1) {
+        MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
+                 tag, comm);
+      }
+      //printf("node %d ordering %d\n",rank,myordering);
+
+      /* receive, reduce, and forward data */
+
+      /* send only */
+      if (myordering == 0) {
+        if (header_buf[myordering + 1] == -1) {
+          to = 0;
+        } else {
+          to = header_buf[myordering + 1];
+        }
+        MPI_Send(rbuf, count, datatype, to, tag, comm);
+      }
+
+      /* recv, reduce, send */
+      else {
+        if (header_buf[myordering + 1] == -1) {
+          to = 0;
+        } else {
+          to = header_buf[myordering + 1];
+        }
+        from = header_buf[myordering - 1];
+        MPI_Recv(tmp_buf, count, datatype, header_buf[myordering - 1], tag,
+                 comm, &status);
+        star_reduction(op, tmp_buf, rbuf, &count, &datatype);
+        MPI_Send(rbuf, count, datatype, to, tag, comm);
+      }
+    }                           /* non-root */
+  }
+  /* pipeline bcast */
+  else {
+    //    printf("node %d start\n",rank);
+
+    send_request_array =
+        (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
+    recv_request_array =
+        (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
+    send_status_array =
+        (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
+    recv_status_array =
+        (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
+
+    if (rank == 0) {
+      sent_count = 0;
+
+      int will_send[MAX_NODE];
+      for (i = 0; i < MAX_NODE; i++)
+        will_send[i] = 0;
+
+      /* loop until all data are received (sent) */
+      while (sent_count < (size - 1)) {
+        int k;
+        for (k = 0; k < 1; k++) {
+          for (i = 1; i < size; i++) {
+            //if (i == rank)
+            //continue;
+            if ((already_received[i] == 0) && (will_send[i] == 0)) {
+              MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
+                         &temp_status_array[i]);
+              if (flag_array[i] == 1) {
+                will_send[i] = 1;
+                MPI_Recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD,
+                         &status);
+                //printf("recv from %d\n",i);
+                i = 1;
+              }
+            }
+          }
+        }                       /* end of probing */
+
+        header_index = 0;
+
+        /* recv 1-byte message */
+        for (i = 1; i < size; i++) {
+          //if (i==rank)
+          //continue;
+          /* message arrived in this round (put in the header) */
+          if ((will_send[i] == 1) && (already_received[i] == 0)) {
+            header_buf[header_index] = i;
+            header_index++;
+            sent_count++;
+
+            /* will send in the next step */
+            already_received[i] = 1;
+          }
+        }
+
+        /* send header followed by data */
+        if (header_index != 0) {
+          header_buf[header_index] = -1;
+          to = header_buf[0];
+
+          /* send header */
+          MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
+
+          /* recv data - pipeline */
+          from = header_buf[header_index - 1];
+          for (i = 0; i < pipe_length; i++) {
+            MPI_Recv(tmp_buf + (i * increment), segment, datatype, from, tag,
+                     comm, &status);
+            star_reduction(op, tmp_buf + (i * increment),
+                           (char *)rbuf + (i * increment), &segment, &datatype);
+          }
+        }
+      }                         /* while loop (sent_count < size-1 ) */
+    }
+
+    /* root */
+    /* none root */
+    else {
+      /* send 1-byte message to root */
+      MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
+
+
+      /* wait for header forward when required */
+      MPI_Irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
+                &request);
+      MPI_Wait(&request, MPI_STATUS_IGNORE);
+
+      /* search for where it is */
+      int myordering = 0;
+
+      while (rank != header_buf[myordering]) {
+        myordering++;
+      }
+
+      /* send header when required */
+      if (header_buf[myordering + 1] != -1) {
+        MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
+                 tag, comm);
+      }
+
+      /* (receive, reduce), and send data */
+      if (header_buf[myordering + 1] == -1) {
+        to = 0;
+      } else {
+        to = header_buf[myordering + 1];
+      }
+
+      /* send only */
+      if (myordering == 0) {
+        for (i = 0; i < pipe_length; i++) {
+          MPI_Isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm,
+                    &send_request_array[i]);
+        }
+        MPI_Waitall((pipe_length), send_request_array, send_status_array);
+      }
+
+      /* receive, reduce, and send */
+      else {
+        from = header_buf[myordering - 1];
+        for (i = 0; i < pipe_length; i++) {
+          MPI_Irecv(tmp_buf + (i * increment), segment, datatype, from, tag,
+                    comm, &recv_request_array[i]);
+        }
+        for (i = 0; i < pipe_length; i++) {
+          MPI_Wait(&recv_request_array[i], MPI_STATUS_IGNORE);
+          star_reduction(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment),
+                         &segment, &datatype);
+          MPI_Isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm,
+                    &send_request_array[i]);
+        }
+        MPI_Waitall((pipe_length), send_request_array, send_status_array);
+      }
+    }                           /* non-root */
+
+
+
+
+    free(send_request_array);
+    free(recv_request_array);
+    free(send_status_array);
+    free(recv_status_array);
+
+    //printf("node %d done\n",rank);
+  }                             /* end pipeline */
+
+
+  /* if root is not zero send root after finished
+     this can be modified to make it faster by using logical src, dst.
+   */
+  if (root != 0) {
+    if (rank == 0) {
+      MPI_Send(rbuf, count, datatype, root, tag, comm);
+    } else if (rank == root) {
+      MPI_Recv(rbuf, count, datatype, 0, tag, comm, &status);
+    }
+  }
+
+
+  /* when count is not divisible by block size, use default BCAST for the remainder */
+  if ((remainder != 0) && (count > segment)) {
+    MPI_Reduce((char *)buf + (pipe_length * increment),
+              (char *)rbuf + (pipe_length * increment), remainder, datatype, op, root,
+               comm);
+  }
+
+  free(tmp_buf);
+
+  return MPI_SUCCESS;
+}