Logo AND Algorithmique Numérique Distribuée

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