Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of framagit.org:simgrid/simgrid
[simgrid.git] / src / smpi / colls / allreduce / allreduce-rab-rdb.cpp
1 /* Copyright (c) 2013-2019. 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.hpp"
8 namespace simgrid{
9 namespace smpi{
10 int Coll_allreduce_rab_rdb::allreduce(void *sbuff, void *rbuff, int count,
11                                       MPI_Datatype dtype, MPI_Op op,
12                                       MPI_Comm comm)
13 {
14   int tag = COLL_TAG_ALLREDUCE;
15   unsigned int mask, pof2, i, recv_idx, last_idx, send_idx, send_cnt;
16   int dst, newrank, rem, newdst,
17       recv_cnt, *cnts, *disps;
18   MPI_Aint extent;
19   MPI_Status status;
20   void *tmp_buf = NULL;
21
22   unsigned int nprocs = comm->size();
23   int rank = comm->rank();
24
25   extent = dtype->get_extent();
26   tmp_buf = (void *) smpi_get_tmp_sendbuffer(count * extent);
27
28   Datatype::copy(sbuff, count, dtype, rbuff, count, dtype);
29
30   // find nearest power-of-two less than or equal to comm_size
31   pof2 = 1;
32   while (pof2 <= nprocs)
33     pof2 <<= 1;
34   pof2 >>= 1;
35
36   rem = nprocs - pof2;
37
38   // In the non-power-of-two case, all even-numbered
39   // processes of rank < 2*rem send their data to
40   // (rank+1). These even-numbered processes no longer
41   // participate in the algorithm until the very end. The
42   // remaining processes form a nice power-of-two.
43
44   if (rank < 2 * rem) {
45     // even
46     if (rank % 2 == 0) {
47
48       Request::send(rbuff, count, dtype, rank + 1, tag, comm);
49
50       // temporarily set the rank to -1 so that this
51       // process does not pariticipate in recursive
52       // doubling
53       newrank = -1;
54     } else                      // odd
55     {
56       Request::recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
57       // do the reduction on received data. since the
58       // ordering is right, it doesn't matter whether
59       // the operation is commutative or not.
60        if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuff, &count, dtype);
61
62       // change the rank
63       newrank = rank / 2;
64     }
65   }
66
67   else                          // rank >= 2 * rem
68     newrank = rank - rem;
69
70   // If op is user-defined or count is less than pof2, use
71   // recursive doubling algorithm. Otherwise do a reduce-scatter
72   // followed by allgather. (If op is user-defined,
73   // derived datatypes are allowed and the user could pass basic
74   // datatypes on one process and derived on another as long as
75   // the type maps are the same. Breaking up derived
76   // datatypes to do the reduce-scatter is tricky, therefore
77   // using recursive doubling in that case.)
78
79   if (newrank != -1) {
80     // do a reduce-scatter followed by allgather. for the
81     // reduce-scatter, calculate the count that each process receives
82     // and the displacement within the buffer
83
84     cnts = (int *) xbt_malloc(pof2 * sizeof(int));
85     disps = (int *) xbt_malloc(pof2 * sizeof(int));
86
87     for (i = 0; i < (pof2 - 1); i++)
88       cnts[i] = count / pof2;
89     cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
90
91     disps[0] = 0;
92     for (i = 1; i < pof2; i++)
93       disps[i] = disps[i - 1] + cnts[i - 1];
94
95     mask = 0x1;
96     send_idx = recv_idx = 0;
97     last_idx = pof2;
98     while (mask < pof2) {
99       newdst = newrank ^ mask;
100       // find real rank of dest
101       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
102
103       send_cnt = recv_cnt = 0;
104       if (newrank < newdst) {
105         send_idx = recv_idx + pof2 / (mask * 2);
106         for (i = send_idx; i < last_idx; i++)
107           send_cnt += cnts[i];
108         for (i = recv_idx; i < send_idx; i++)
109           recv_cnt += cnts[i];
110       } else {
111         recv_idx = send_idx + pof2 / (mask * 2);
112         for (i = send_idx; i < recv_idx; i++)
113           send_cnt += cnts[i];
114         for (i = recv_idx; i < last_idx; i++)
115           recv_cnt += cnts[i];
116       }
117
118       // Send data from recvbuf. Recv into tmp_buf
119       Request::sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
120                    dtype, dst, tag,
121                    (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
122                    dtype, dst, tag, comm, &status);
123
124       // tmp_buf contains data received in this step.
125       // recvbuf contains data accumulated so far
126
127       // This algorithm is used only for predefined ops
128       // and predefined ops are always commutative.
129       if(op!=MPI_OP_NULL) op->apply( (char *) tmp_buf + disps[recv_idx] * extent,
130                         (char *) rbuff + disps[recv_idx] * extent, &recv_cnt, dtype);
131
132       // update send_idx for next iteration
133       send_idx = recv_idx;
134       mask <<= 1;
135
136       // update last_idx, but not in last iteration because the value
137       // is needed in the allgather step below.
138       if (mask < pof2)
139         last_idx = recv_idx + pof2 / mask;
140     }
141
142     // now do the allgather
143
144     mask >>= 1;
145     while (mask > 0) {
146       newdst = newrank ^ mask;
147       // find real rank of dest
148       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
149
150       send_cnt = recv_cnt = 0;
151       if (newrank < newdst) {
152         // update last_idx except on first iteration
153         if (mask != pof2 / 2)
154           last_idx = last_idx + pof2 / (mask * 2);
155
156         recv_idx = send_idx + pof2 / (mask * 2);
157         for (i = send_idx; i < recv_idx; i++)
158           send_cnt += cnts[i];
159         for (i = recv_idx; i < last_idx; i++)
160           recv_cnt += cnts[i];
161       } else {
162         recv_idx = send_idx - pof2 / (mask * 2);
163         for (i = send_idx; i < last_idx; i++)
164           send_cnt += cnts[i];
165         for (i = recv_idx; i < send_idx; i++)
166           recv_cnt += cnts[i];
167       }
168
169       Request::sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
170                    dtype, dst, tag,
171                    (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
172                    dtype, dst, tag, comm, &status);
173
174       if (newrank > newdst)
175         send_idx = recv_idx;
176
177       mask >>= 1;
178     }
179
180     free(cnts);
181     free(disps);
182
183   }
184   // In the non-power-of-two case, all odd-numbered processes of
185   // rank < 2 * rem send the result to (rank-1), the ranks who didn't
186   // participate above.
187
188   if (rank < 2 * rem) {
189     if (rank % 2)               // odd
190       Request::send(rbuff, count, dtype, rank - 1, tag, comm);
191     else                        // even
192       Request::recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
193   }
194
195   smpi_free_tmp_buffer(tmp_buf);
196   return MPI_SUCCESS;
197 }
198 }
199 }