Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[DOC] Fixed even more errors.
[simgrid.git] / src / smpi / colls / alltoall-ring-light-barrier.c
index f901445..eb3de84 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"
 /*****************************************************************************
 
  * Function: alltoall_ring_light_barrier
@@ -31,21 +37,21 @@ smpi_coll_tuned_alltoall_ring_light_barrier(void *send_buff, int send_count,
   MPI_Aint send_chunk, recv_chunk;
   MPI_Status s;
   int i, src, dst, rank, num_procs, next_dst, next_src;
-  int tag = 1, success = 1;     /*, failure = 0; */
+  int tag = COLL_TAG_ALLTOALL;
 
   char send_sync = 'a', recv_sync = 'b';
   char *send_ptr = (char *) send_buff;
   char *recv_ptr = (char *) recv_buff;
 
-  MPI_Comm_rank(comm, &rank);
-  MPI_Comm_size(comm, &num_procs);
-  MPI_Type_extent(send_type, &send_chunk);
-  MPI_Type_extent(recv_type, &recv_chunk);
+  rank = smpi_comm_rank(comm);
+  num_procs = smpi_comm_size(comm);
+  send_chunk = smpi_datatype_get_extent(send_type);
+  recv_chunk = smpi_datatype_get_extent(recv_type);
 
   send_chunk *= send_count;
   recv_chunk *= recv_count;
 
-  MPI_Sendrecv(send_ptr + rank * send_chunk, send_count, send_type, rank, tag,
+  smpi_mpi_sendrecv(send_ptr + rank * send_chunk, send_count, send_type, rank, tag,
                recv_ptr + rank * recv_chunk, recv_count, recv_type, rank, tag,
                comm, &s);
 
@@ -53,18 +59,18 @@ smpi_coll_tuned_alltoall_ring_light_barrier(void *send_buff, int send_count,
     src = (rank - i + num_procs) % num_procs;
     dst = (rank + i) % num_procs;
 
-    MPI_Sendrecv(send_ptr + dst * send_chunk, send_count, send_type,
+    smpi_mpi_sendrecv(send_ptr + dst * send_chunk, send_count, send_type,
                  dst, tag, recv_ptr + src * recv_chunk, recv_count,
                  recv_type, src, tag, comm, &s);
 
     if ((i + 1) < num_procs) {
       next_src = (rank - (i + 1) + num_procs) % num_procs;
       next_dst = (rank + (i + 1) + num_procs) % num_procs;
-      MPI_Sendrecv(&send_sync, 1, MPI_CHAR, next_src, tag,
+      smpi_mpi_sendrecv(&send_sync, 1, MPI_CHAR, next_src, tag,
                    &recv_sync, 1, MPI_CHAR, next_dst, tag, comm, &s);
 
     }
   }
 
-  return success;
+  return MPI_SUCCESS;
 }