Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Use simgrid function instead of MPI in collectives
[simgrid.git] / src / smpi / colls / bcast-arrival-nb.c
index a4246c4..9ff27b4 100644 (file)
@@ -1,4 +1,4 @@
-#include "colls.h"
+#include "colls_private.h"
 
 static int bcast_NTSL_segment_size_in_byte = 8192;
 
@@ -32,15 +32,15 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
   char temp_buf[MAX_NODE];
 
   MPI_Aint extent;
-  MPI_Type_extent(datatype, &extent);
+  extent = smpi_datatype_get_extent(datatype);
 
   /* destination */
   int to;
 
 
 
-  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
-  MPI_Comm_size(MPI_COMM_WORLD, &size);
+  rank = smpi_comm_rank(MPI_COMM_WORLD);
+  size = smpi_comm_size(MPI_COMM_WORLD);
 
 
   /* segment is segment size in number of elements (not bytes) */
@@ -61,9 +61,9 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
    */
   if (root != 0) {
     if (rank == root) {
-      MPI_Send(buf, count, datatype, 0, tag, comm);
+      smpi_mpi_send(buf, count, datatype, 0, tag, comm);
     } else if (rank == 0) {
-      MPI_Recv(buf, count, datatype, root, tag, comm, &status);
+      smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
     }
   }
 
@@ -83,7 +83,7 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
         //      for (j=0;j<1000;j++) {
         for (i = 1; i < size; i++) {
           if (already_sent[i] == 0)
-            MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
+            smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
                        MPI_STATUSES_IGNORE);
         }
         //}
@@ -94,7 +94,7 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
           /* message arrive */
           if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
-            MPI_Recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
+            smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
             header_buf[header_index] = i;
             header_index++;
             sent_count++;
@@ -108,8 +108,8 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
         if (header_index != 0) {
           header_buf[header_index] = -1;
           to = header_buf[0];
-          MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
-          MPI_Send(buf, count, datatype, to, tag, comm);
+          smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
+          smpi_mpi_send(buf, count, datatype, to, tag, comm);
         }
 
         /* randomly MPI_Send to one */
@@ -119,8 +119,8 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
             if (already_sent[i] == 0) {
               header_buf[0] = i;
               header_buf[1] = -1;
-              MPI_Send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
-              MPI_Send(buf, count, datatype, i, tag, comm);
+              smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
+              smpi_mpi_send(buf, count, datatype, i, tag, comm);
               already_sent[i] = 1;
               sent_count++;
               break;
@@ -136,12 +136,12 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
     else {
 
       /* send 1-byte message to root */
-      MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
+      smpi_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,
+      smpi_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);
+      smpi_mpi_recv(buf, count, datatype, MPI_ANY_SOURCE, tag, comm, &status);
 
       /* search for where it is */
       int myordering = 0;
@@ -151,22 +151,22 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
       /* send header followed by data */
       if (header_buf[myordering + 1] != -1) {
-        MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
+        smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
                  tag, comm);
-        MPI_Send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
+        smpi_mpi_send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
       }
     }
   }
   /* pipeline bcast */
   else {
     send_request_array =
-        (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
+        (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
     recv_request_array =
-        (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
+        (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
     send_status_array =
-        (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
+        (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
     recv_status_array =
-        (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
+        (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
 
     if (rank == 0) {
       sent_count = 0;
@@ -183,11 +183,11 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
         for (k = 0; k < 3; k++) {
           for (i = 1; i < size; i++) {
             if ((already_sent[i] == 0) && (will_send[i] == 0)) {
-              MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
+              smpi_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,
+                smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD,
                          &status);
                 i = 1;
               }
@@ -238,7 +238,7 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
           //start = MPI_Wtime();
 
           /* send header */
-          MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
+          smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
 
           //total = MPI_Wtime() - start;
           //total *= 1000;
@@ -250,16 +250,16 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
           if (0 == 1) {
             //if (header_index == 1) {
-            MPI_Send(buf, count, datatype, to, tag, comm);
+            smpi_mpi_send(buf, count, datatype, to, tag, comm);
           }
 
 
           /* send data - pipeline */
           else {
             for (i = 0; i < pipe_length; i++) {
-              MPI_Send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
+              smpi_mpi_send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
             }
-            //MPI_Waitall((pipe_length), send_request_array, send_status_array);
+            //smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
           }
           //total = MPI_Wtime() - start;
           //total *= 1000;
@@ -279,16 +279,16 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
               to = i;
 
               //start = MPI_Wtime();
-              MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
+              smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
 
               /* still need to chop data so that we can use the same non-root code */
               for (j = 0; j < pipe_length; j++) {
-                MPI_Send((char *)buf + (j * increment), segment, datatype, to, tag,
+                smpi_mpi_send((char *)buf + (j * increment), segment, datatype, to, tag,
                          comm);
               }
 
-              //MPI_Send(buf,count,datatype,to,tag,comm);
-              //MPI_Wait(&request,MPI_STATUS_IGNORE);
+              //smpi_mpi_send(buf,count,datatype,to,tag,comm);
+              //smpi_mpi_wait(&request,MPI_STATUS_IGNORE);
 
               //total = MPI_Wtime() - start;
               //total *= 1000;
@@ -310,10 +310,10 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
       /* probe before exit in case there are messages to recv */
       for (i = 1; i < size; i++) {
-        MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
+        smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
                    &temp_status_array[i]);
         if (flag_array[i] == 1)
-          MPI_Recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
+          smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
       }
     }
 
@@ -322,16 +322,15 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
     else {
 
       /* if root already send a message to this node, don't send one-byte message */
-      MPI_Iprobe(0, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[0], &status);
+      smpi_mpi_iprobe(0, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[0], &status);
 
       /* send 1-byte message to root */
       if (flag_array[0] == 0)
-        MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
+        smpi_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);
+      request = smpi_mpi_irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm);
+      smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
 
       /* search for where it is */
       int myordering = 0;
@@ -341,7 +340,7 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
       /* send header when required */
       if (header_buf[myordering + 1] != -1) {
-        MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
+        smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
                  tag, comm);
       }
 
@@ -349,25 +348,24 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
       if (0 == -1) {
         //if (header_buf[1] == -1) {
-        MPI_Irecv(buf, count, datatype, 0, tag, comm, &request);
-        MPI_Wait(&request, MPI_STATUS_IGNORE);
+        request = smpi_mpi_irecv(buf, count, datatype, 0, tag, comm);
+        smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
         //printf("\t\tnode %d ordering = %d receive data from root\n",rank,myordering);
       } else {
         for (i = 0; i < pipe_length; i++) {
-          MPI_Irecv((char *)buf + (i * increment), segment, datatype, MPI_ANY_SOURCE,
-                    tag, comm, &recv_request_array[i]);
+          recv_request_array[i] = smpi_mpi_irecv((char *)buf + (i * increment), segment, datatype, MPI_ANY_SOURCE,
+                    tag, comm);
         }
       }
 
       /* send data */
       if (header_buf[myordering + 1] != -1) {
         for (i = 0; i < pipe_length; i++) {
-          MPI_Wait(&recv_request_array[i], MPI_STATUS_IGNORE);
-          MPI_Isend((char *)buf + (i * increment), segment, datatype,
-                    header_buf[myordering + 1], tag, comm,
-                    &send_request_array[i]);
+          smpi_mpi_wait(&recv_request_array[i], MPI_STATUS_IGNORE);
+          send_request_array[i] = smpi_mpi_isend((char *)buf + (i * increment), segment, datatype,
+                    header_buf[myordering + 1], tag, comm);
         }
-        MPI_Waitall((pipe_length), send_request_array, send_status_array);
+        smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
       }
 
     }
@@ -380,7 +378,8 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
   /* when count is not divisible by block size, use default BCAST for the remainder */
   if ((remainder != 0) && (count > segment)) {
-    MPI_Bcast((char *)buf + (pipe_length * increment), remainder, datatype, root, comm);
+    XBT_WARN("MPI_bcast_arrival_nb use default MPI_bcast.");     
+    smpi_mpi_bcast((char *)buf + (pipe_length * increment), remainder, datatype, root, comm);
   }
 
   return MPI_SUCCESS;