Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add tesh files to test all new collectives
[simgrid.git] / src / smpi / colls / alltoall-2dmesh.c
index 2b1a2a0..644e6b9 100644 (file)
@@ -27,7 +27,7 @@ XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_colls, smpi,
  * Auther: Ahmad Faraj
 
 ****************************************************************************/
-int alltoall_check_is_2dmesh(int num, int *i, int *j)
+static int alltoall_check_is_2dmesh(int num, int *i, int *j)
 {
   int x, max = num / 2;
   x = sqrt(num);
@@ -53,8 +53,7 @@ int alltoall_check_is_2dmesh(int num, int *i, int *j)
 int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
                                     MPI_Datatype send_type,
                                     void *recv_buff, int recv_count,
-                                    MPI_Datatype recv_type,
-                                    MPI_Comm comm)
+                                    MPI_Datatype recv_type, MPI_Comm comm)
 {
   MPI_Status *statuses, s;
   MPI_Request *reqs, *req_ptr;;
@@ -62,18 +61,16 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
 
   char *tmp_buff1, *tmp_buff2;
   int i, j, src, dst, rank, num_procs, count, num_reqs;
-  int rows, cols, my_row, my_col, X, Y, send_offset, recv_offset;
-  int two_dsize, my_row_base, my_col_base, src_row_base, block_size;
-  int tag = 1, failure = 0, success = 1;
+  int X, Y, send_offset, recv_offset;
+  int my_row_base, my_col_base, src_row_base, block_size;
+  int tag = 1;
 
   MPI_Comm_rank(comm, &rank);
   MPI_Comm_size(comm, &num_procs);
   MPI_Type_extent(send_type, &extent);
 
   if (!alltoall_check_is_2dmesh(num_procs, &X, &Y))
-    return failure;
-
-  two_dsize = X * Y;
+    return MPI_ERR_OTHER;
 
   my_row_base = (rank / Y) * Y;
   my_col_base = rank % Y;
@@ -84,14 +81,14 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
   if (!tmp_buff1) {
     XBT_DEBUG("alltoall-2dmesh_shoot.c:88: cannot allocate memory");
     MPI_Finalize();
-    exit(failure);
+    exit(MPI_ERR_OTHER);
   }
 
   tmp_buff2 = (char *) malloc(block_size * Y);
   if (!tmp_buff2) {
     XBT_WARN("alltoall-2dmesh_shoot.c:88: cannot allocate memory");
     MPI_Finalize();
-    exit(failure);
+    exit(MPI_ERR_OTHER);
   }
 
 
@@ -105,7 +102,7 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
   if (!reqs) {
     XBT_WARN("alltoall-2dmesh_shoot.c:88: cannot allocate memory");
     MPI_Finalize();
-    exit(failure);
+    exit(MPI_ERR_OTHER);
   }
 
   req_ptr = reqs;
@@ -139,14 +136,15 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
     recv_offset = (my_row_base * block_size) + (i * block_size);
 
     if (i + my_row_base == rank)
-      MPI_Sendrecv(send_buff + recv_offset, send_count, send_type,
-                   rank, tag, recv_buff + recv_offset, recv_count,
-                   recv_type, rank, tag, comm, &s);
+      MPI_Sendrecv((char *) send_buff + recv_offset, send_count, send_type,
+                   rank, tag,
+                   (char *) recv_buff + recv_offset, recv_count, recv_type,
+                   rank, tag, comm, &s);
 
     else
       MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type,
                    rank, tag,
-                   recv_buff + recv_offset, recv_count, recv_type,
+                   (char *) recv_buff + recv_offset, recv_count, recv_type,
                    rank, tag, comm, &s);
   }
 
@@ -157,7 +155,7 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
       continue;
     src_row_base = (src / Y) * Y;
 
-    MPI_Irecv(recv_buff + src_row_base * block_size, recv_count * Y,
+    MPI_Irecv((char *) recv_buff + src_row_base * block_size, recv_count * Y,
               recv_type, src, tag, comm, req_ptr++);
   }
 
@@ -171,10 +169,9 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
       send_offset = (dst + j * num_procs) * block_size;
 
       if (j + my_row_base == rank)
-        MPI_Sendrecv(send_buff + dst * block_size, send_count, send_type,
-                     rank, tag,
-                     tmp_buff2 + recv_offset, recv_count, recv_type,
-                     rank, tag, comm, &s);
+        MPI_Sendrecv((char *) send_buff + dst * block_size, send_count,
+                     send_type, rank, tag, tmp_buff2 + recv_offset, recv_count,
+                     recv_type, rank, tag, comm, &s);
       else
         MPI_Sendrecv(tmp_buff1 + send_offset, send_count, send_type,
                      rank, tag,
@@ -191,5 +188,5 @@ int smpi_coll_tuned_alltoall_2dmesh(void *send_buff, int send_count,
   free(statuses);
   free(tmp_buff1);
   free(tmp_buff2);
-  return success;
+  return MPI_SUCCESS;
 }