Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Use simgrid function instead of MPI in collectives
[simgrid.git] / src / smpi / colls / alltoall-3dmesh.c
index 7ffac2b..ca10de7 100644 (file)
@@ -1,4 +1,4 @@
-#include "colls.h"
+#include "colls_private.h"
 #include <math.h>
 
 /*****************************************************************************
@@ -54,9 +54,9 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
 
   char *tmp_buff1, *tmp_buff2;
 
-  MPI_Comm_rank(comm, &rank);
-  MPI_Comm_size(comm, &num_procs);
-  MPI_Type_extent(send_type, &extent);
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  extent = smpi_datatype_get_extent(send_type);
 
   if (!alltoall_check_is_3dmesh(num_procs, &X, &Y, &Z))
     return MPI_ERR_OTHER;
@@ -86,7 +86,7 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
 
   send_offset = recv_offset = (rank % two_dsize) * block_size * num_procs;
 
-  MPI_Sendrecv(send_buff, send_count * num_procs, send_type, rank, tag,
+  smpi_mpi_sendrecv(send_buff, send_count * num_procs, send_type, rank, tag,
                tmp_buff1 + recv_offset, num_procs * recv_count,
                recv_type, rank, tag, comm, &status);
 
@@ -97,18 +97,17 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
     if (src == rank)
       continue;
     recv_offset = (src % two_dsize) * block_size * num_procs;
-    MPI_Irecv(tmp_buff1 + recv_offset, count, recv_type, src, tag, comm,
-              req_ptr++);
+    *(req_ptr++) = smpi_mpi_irecv(tmp_buff1 + recv_offset, count, recv_type, src, tag, comm);
   }
 
   for (i = 0; i < Y; i++) {
     dst = i + my_row_base;
     if (dst == rank)
       continue;
-    MPI_Send(send_buff, count, send_type, dst, tag, comm);
+    smpi_mpi_send(send_buff, count, send_type, dst, tag, comm);
   }
 
-  MPI_Waitall(Y - 1, reqs, statuses);
+  smpi_mpi_waitall(Y - 1, reqs, statuses);
   req_ptr = reqs;
 
 
@@ -120,8 +119,8 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
     src_row_base = (src / X) * X;
 
     recv_offset = (src_row_base % two_dsize) * block_size * num_procs;
-    MPI_Irecv(tmp_buff1 + recv_offset, recv_count * num_procs * Y,
-              recv_type, src, tag, comm, req_ptr++);
+    *(req_ptr++) = smpi_mpi_irecv(tmp_buff1 + recv_offset, recv_count * num_procs * Y,
+              recv_type, src, tag, comm);
   }
 
   send_offset = (my_row_base % two_dsize) * block_size * num_procs;
@@ -129,17 +128,17 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
     dst = (i * Y + my_col_base);
     if (dst == rank)
       continue;
-    MPI_Send(tmp_buff1 + send_offset, send_count * num_procs * Y, send_type,
+    smpi_mpi_send(tmp_buff1 + send_offset, send_count * num_procs * Y, send_type,
              dst, tag, comm);
   }
 
-  MPI_Waitall(X - 1, reqs, statuses);
+  smpi_mpi_waitall(X - 1, reqs, statuses);
   req_ptr = reqs;
 
   for (i = 0; i < two_dsize; i++) {
     send_offset = (rank * block_size) + (i * block_size * num_procs);
     recv_offset = (my_z_base * block_size) + (i * block_size);
-    MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type, rank, tag,
+    smpi_mpi_sendrecv(tmp_buff1 + send_offset, send_count, send_type, rank, tag,
                  (char *) recv_buff + recv_offset, recv_count, recv_type,
                  rank, tag, comm, &status);
   }
@@ -150,8 +149,8 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
 
     recv_offset = (src_z_base * block_size);
 
-    MPI_Irecv((char *) recv_buff + recv_offset, recv_count * two_dsize,
-              recv_type, src, tag, comm, req_ptr++);
+    *(req_ptr++) = smpi_mpi_irecv((char *) recv_buff + recv_offset, recv_count * two_dsize,
+              recv_type, src, tag, comm);
   }
 
   for (i = 1; i < Z; i++) {
@@ -160,18 +159,18 @@ int smpi_coll_tuned_alltoall_3dmesh(void *send_buff, int send_count,
     recv_offset = 0;
     for (j = 0; j < two_dsize; j++) {
       send_offset = (dst + j * num_procs) * block_size;
-      MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type,
+      smpi_mpi_sendrecv(tmp_buff1 + send_offset, send_count, send_type,
                    rank, tag, tmp_buff2 + recv_offset, recv_count,
                    recv_type, rank, tag, comm, &status);
 
       recv_offset += block_size;
     }
 
-    MPI_Send(tmp_buff2, send_count * two_dsize, send_type, dst, tag, comm);
+    smpi_mpi_send(tmp_buff2, send_count * two_dsize, send_type, dst, tag, comm);
 
   }
 
-  MPI_Waitall(Z - 1, reqs, statuses);
+  smpi_mpi_waitall(Z - 1, reqs, statuses);
 
   free(reqs);
   free(statuses);