Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove warning
[simgrid.git] / src / smpi / colls / reduce-binomial.c
1 #include "colls_private.h"
2
3 //#include <star-reduction.c>
4
5 int smpi_coll_tuned_reduce_binomial(void *sendbuf, void *recvbuf, int count,
6                                     MPI_Datatype datatype, MPI_Op op, int root,
7                                     MPI_Comm comm)
8 {
9   MPI_Status status;
10   int comm_size, rank;
11   int mask, relrank, source;
12   int dst;
13   int tag = 4321;
14   MPI_Aint extent;
15   void *tmp_buf;
16
17   if (count == 0)
18     return 0;
19   rank = smpi_comm_rank(comm);
20   comm_size = smpi_comm_size(comm);
21
22   extent = smpi_datatype_get_extent(datatype);
23
24   tmp_buf = (void *) xbt_malloc(count * extent);
25
26   smpi_mpi_sendrecv(sendbuf, count, datatype, rank, tag,
27                recvbuf, count, datatype, rank, tag, comm, &status);
28   mask = 1;
29   relrank = (rank - root + comm_size) % comm_size;
30
31   while (mask < comm_size) {
32     /* Receive */
33     if ((mask & relrank) == 0) {
34       source = (relrank | mask);
35       if (source < comm_size) {
36         source = (source + root) % comm_size;
37         smpi_mpi_recv(tmp_buf, count, datatype, source, tag, comm, &status);
38         smpi_op_apply(op, tmp_buf, recvbuf, &count, &datatype);
39       }
40     } else {
41       dst = ((relrank & (~mask)) + root) % comm_size;
42       smpi_mpi_send(recvbuf, count, datatype, dst, tag, comm);
43       break;
44     }
45     mask <<= 1;
46   }
47
48   free(tmp_buf);
49
50   return 0;
51 }