Logo AND Algorithmique Numérique Distribuée

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