Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'smpi_cpp'
[simgrid.git] / src / smpi / colls / allreduce-rdb.cpp
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "colls_private.h"
8 //#include <star-reduction.c>
9
10 int smpi_coll_tuned_allreduce_rdb(void *sbuff, void *rbuff, int count,
11                                   MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)
12 {
13   int nprocs, rank, tag = COLL_TAG_ALLREDUCE;
14   int mask, dst, pof2, newrank, rem, newdst;
15   MPI_Aint extent, lb;
16   MPI_Status status;
17   void *tmp_buf = NULL;
18   /*
19      #ifdef MPICH2_REDUCTION
20      MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
21      #else
22      MPI_User_function *uop;
23      struct MPIR_OP *op_ptr;
24      op_ptr = MPIR_ToPointer(op);
25      uop  = op_ptr->op;
26      #endif
27    */
28   nprocs=comm->size();
29   rank=comm->rank();
30
31   smpi_datatype_extent(dtype, &lb, &extent);
32   tmp_buf = (void *) smpi_get_tmp_sendbuffer(count * extent);
33
34   smpi_mpi_sendrecv(sbuff, count, dtype, rank, 500,
35                rbuff, count, dtype, rank, 500, comm, &status);
36
37   // find nearest power-of-two less than or equal to comm_size
38   pof2 = 1;
39   while (pof2 <= nprocs)
40     pof2 <<= 1;
41   pof2 >>= 1;
42
43   rem = nprocs - pof2;
44
45   // In the non-power-of-two case, all even-numbered
46   // processes of rank < 2*rem send their data to
47   // (rank+1). These even-numbered processes no longer
48   // participate in the algorithm until the very end. The
49   // remaining processes form a nice power-of-two. 
50
51   if (rank < 2 * rem) {
52     // even       
53     if (rank % 2 == 0) {
54
55       smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
56
57       // temporarily set the rank to -1 so that this
58       // process does not pariticipate in recursive
59       // doubling
60       newrank = -1;
61     } else                      // odd
62     {
63       smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
64       // do the reduction on received data. since the
65       // ordering is right, it doesn't matter whether
66       // the operation is commutative or not.
67       smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
68
69       // change the rank 
70       newrank = rank / 2;
71     }
72   }
73
74   else                          // rank >= 2 * rem 
75     newrank = rank - rem;
76
77   // If op is user-defined or count is less than pof2, use
78   // recursive doubling algorithm. Otherwise do a reduce-scatter
79   // followed by allgather. (If op is user-defined,
80   // derived datatypes are allowed and the user could pass basic
81   // datatypes on one process and derived on another as long as
82   // the type maps are the same. Breaking up derived
83   // datatypes to do the reduce-scatter is tricky, therefore
84   // using recursive doubling in that case.) 
85
86   if (newrank != -1) {
87     mask = 0x1;
88     while (mask < pof2) {
89       newdst = newrank ^ mask;
90       // find real rank of dest 
91       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
92
93       // Send the most current data, which is in recvbuf. Recv
94       // into tmp_buf 
95       smpi_mpi_sendrecv(rbuff, count, dtype, dst, tag, tmp_buf, count, dtype,
96                    dst, tag, comm, &status);
97
98       // tmp_buf contains data received in this step.
99       // recvbuf contains data accumulated so far 
100
101       // op is commutative OR the order is already right
102       // we assume it is commuttive op
103       //      if (op -> op_commute  || (dst < rank))
104       if ((dst < rank)) {
105         smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
106       } else                    // op is noncommutative and the order is not right
107       {
108         smpi_op_apply(op, rbuff, tmp_buf, &count, &dtype);
109
110         // copy result back into recvbuf
111         smpi_mpi_sendrecv(tmp_buf, count, dtype, rank, tag, rbuff, count,
112                      dtype, rank, tag, comm, &status);
113       }
114       mask <<= 1;
115     }
116   }
117   // In the non-power-of-two case, all odd-numbered processes of
118   // rank < 2 * rem send the result to (rank-1), the ranks who didn't
119   // participate above.
120
121   if (rank < 2 * rem) {
122     if (rank % 2)               // odd 
123       smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
124     else                        // even 
125       smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
126   }
127
128   smpi_free_tmp_buffer(tmp_buf);
129   return MPI_SUCCESS;
130 }