Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright notices
[simgrid.git] / src / smpi / colls / allgather-3dmesh.c
index 035e981..ab964a1 100644 (file)
@@ -1,4 +1,10 @@
-#include "colls.h"
+/* Copyright (c) 2013-2014. The SimGrid Team.
+ * All rights reserved.                                                     */
+
+/* This program is free software; you can redistribute it and/or modify it
+ * under the terms of the license (GNU LGPL) which comes with this package. */
+
+#include "colls_private.h"
 
 /*****************************************************************************
 
@@ -97,15 +103,15 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
   int i, src, dst, rank, num_procs, block_size, my_z_base;
   int my_z, X, Y, Z, send_offset, recv_offset;
   int two_dsize, my_row_base, my_col_base, src_row_base, src_z_base, num_reqs;
-  int success = 0;
-  int failure = 1;
-  int tag = 1;
+  int tag = COLL_TAG_ALLGATHER;
 
-  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 (!is_3dmesh(num_procs, &X, &Y, &Z))
+    THROWF(arg_error,0, "allgather_3dmesh algorithm can't be used with this number of processes! ");
 
-  is_3dmesh(num_procs, &X, &Y, &Z);
 
   num_reqs = X;
 
@@ -123,18 +129,13 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
 
   block_size = extent * send_count;
 
-  req = (MPI_Request *) malloc(num_reqs * sizeof(MPI_Request));
-  if (!req) {
-    printf("allgather-3dmesh-shoot.c:85: cannot allocate memory\n");
-    MPI_Finalize();
-    exit(failure);
-  }
+  req = (MPI_Request *) xbt_malloc(num_reqs * sizeof(MPI_Request));
 
   req_ptr = req;
 
   // do local allgather/local copy 
   recv_offset = rank * block_size;
-  MPIR_Localcopy(send_buff, send_count, send_type, (char *)recv_buff + recv_offset,
+  smpi_datatype_copy(send_buff, send_count, send_type, (char *)recv_buff + recv_offset,
                  recv_count, recv_type);
 
   // do rowwise comm 
@@ -143,18 +144,18 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
     if (src == rank)
       continue;
     recv_offset = src * block_size;
-    MPIC_Irecv((char *)recv_buff + recv_offset, send_count, recv_type, src, tag,
-               comm, req_ptr++);
+    *(req_ptr++) = smpi_mpi_irecv((char *)recv_buff + recv_offset, send_count, recv_type, src, tag,
+               comm);
   }
 
   for (i = 0; i < Y; i++) {
     dst = i + my_row_base;
     if (dst == rank)
       continue;
-    MPIC_Send(send_buff, send_count, send_type, dst, tag, comm);
+    smpi_mpi_send(send_buff, send_count, send_type, dst, tag, comm);
   }
 
-  MPI_Waitall(Y - 1, req, MPI_STATUSES_IGNORE);
+  smpi_mpi_waitall(Y - 1, req, MPI_STATUSES_IGNORE);
   req_ptr = req;
 
   // do colwise comm, it does not matter here if i*X or i *Y since X == Y
@@ -166,8 +167,8 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
 
     src_row_base = (src / X) * X;
     recv_offset = src_row_base * block_size;
-    MPIC_Irecv((char *)recv_buff + recv_offset, recv_count * Y, recv_type, src, tag,
-               comm, req_ptr++);
+    *(req_ptr++) = smpi_mpi_irecv((char *)recv_buff + recv_offset, recv_count * Y, recv_type, src, tag,
+               comm);
   }
 
   send_offset = my_row_base * block_size;
@@ -176,11 +177,11 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
     dst = (i * Y + my_col_base);
     if (dst == rank)
       continue;
-    MPIC_Send((char *)recv_buff + send_offset, send_count * Y, send_type, dst, tag,
+    smpi_mpi_send((char *)recv_buff + send_offset, send_count * Y, send_type, dst, tag,
               comm);
   }
 
-  MPI_Waitall(X - 1, req, MPI_STATUSES_IGNORE);
+  smpi_mpi_waitall(X - 1, req, MPI_STATUSES_IGNORE);
   req_ptr = req;
 
   for (i = 1; i < Z; i++) {
@@ -189,19 +190,19 @@ int smpi_coll_tuned_allgather_3dmesh(void *send_buff, int send_count,
 
     recv_offset = (src_z_base * block_size);
 
-    MPIC_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++) {
     dst = (rank + i * two_dsize) % num_procs;
     send_offset = my_z_base * block_size;
-    MPIC_Send((char *)recv_buff + send_offset, send_count * two_dsize, send_type,
+    smpi_mpi_send((char *)recv_buff + send_offset, send_count * two_dsize, send_type,
               dst, tag, comm);
   }
-  MPI_Waitall(Z - 1, req, MPI_STATUSES_IGNORE);
+  smpi_mpi_waitall(Z - 1, req, MPI_STATUSES_IGNORE);
 
   free(req);
 
-  return success;
+  return MPI_SUCCESS;
 }