1 /* Copyright (c) 2013-2019. The SimGrid Team.
2 * All rights reserved. */
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. */
7 #include "../colls_private.hpp"
10 int Coll_allreduce_rab_rdb::allreduce(const void *sbuff, void *rbuff, int count,
11 MPI_Datatype dtype, MPI_Op op,
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, recv_cnt;
20 unsigned int nprocs = comm->size();
21 int rank = comm->rank();
23 extent = dtype->get_extent();
24 unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
26 Datatype::copy(sbuff, count, dtype, rbuff, count, dtype);
28 // find nearest power-of-two less than or equal to comm_size
30 while (pof2 <= nprocs)
36 // In the non-power-of-two case, all even-numbered
37 // processes of rank < 2*rem send their data to
38 // (rank+1). These even-numbered processes no longer
39 // participate in the algorithm until the very end. The
40 // remaining processes form a nice power-of-two.
46 Request::send(rbuff, count, dtype, rank + 1, tag, comm);
48 // temporarily set the rank to -1 so that this
49 // process does not pariticipate in recursive
54 Request::recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
55 // do the reduction on received data. since the
56 // ordering is right, it doesn't matter whether
57 // the operation is commutative or not.
58 if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuff, &count, dtype);
65 else // rank >= 2 * rem
68 // If op is user-defined or count is less than pof2, use
69 // recursive doubling algorithm. Otherwise do a reduce-scatter
70 // followed by allgather. (If op is user-defined,
71 // derived datatypes are allowed and the user could pass basic
72 // datatypes on one process and derived on another as long as
73 // the type maps are the same. Breaking up derived
74 // datatypes to do the reduce-scatter is tricky, therefore
75 // using recursive doubling in that case.)
78 // do a reduce-scatter followed by allgather. for the
79 // reduce-scatter, calculate the count that each process receives
80 // and the displacement within the buffer
82 int* cnts = new int[pof2];
83 int* disps = new int[pof2];
85 for (i = 0; i < (pof2 - 1); i++)
86 cnts[i] = count / pof2;
87 cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
90 for (i = 1; i < pof2; i++)
91 disps[i] = disps[i - 1] + cnts[i - 1];
94 send_idx = recv_idx = 0;
97 newdst = newrank ^ mask;
98 // find real rank of dest
99 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
101 send_cnt = recv_cnt = 0;
102 if (newrank < newdst) {
103 send_idx = recv_idx + pof2 / (mask * 2);
104 for (i = send_idx; i < last_idx; i++)
106 for (i = recv_idx; i < send_idx; i++)
109 recv_idx = send_idx + pof2 / (mask * 2);
110 for (i = send_idx; i < recv_idx; i++)
112 for (i = recv_idx; i < last_idx; i++)
116 // Send data from recvbuf. Recv into tmp_buf
117 Request::sendrecv(static_cast<char*>(rbuff) + disps[send_idx] * extent, send_cnt, dtype, dst, tag,
118 tmp_buf + disps[recv_idx] * extent, recv_cnt, dtype, dst, tag, comm, &status);
120 // tmp_buf contains data received in this step.
121 // recvbuf contains data accumulated so far
123 // This algorithm is used only for predefined ops
124 // and predefined ops are always commutative.
125 if (op != MPI_OP_NULL)
126 op->apply(tmp_buf + disps[recv_idx] * extent, static_cast<char*>(rbuff) + disps[recv_idx] * extent, &recv_cnt,
129 // update send_idx for next iteration
133 // update last_idx, but not in last iteration because the value
134 // is needed in the allgather step below.
136 last_idx = recv_idx + pof2 / mask;
139 // now do the allgather
143 newdst = newrank ^ mask;
144 // find real rank of dest
145 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
147 send_cnt = recv_cnt = 0;
148 if (newrank < newdst) {
149 // update last_idx except on first iteration
150 if (mask != pof2 / 2)
151 last_idx = last_idx + pof2 / (mask * 2);
153 recv_idx = send_idx + pof2 / (mask * 2);
154 for (i = send_idx; i < recv_idx; i++)
156 for (i = recv_idx; i < last_idx; i++)
159 recv_idx = send_idx - pof2 / (mask * 2);
160 for (i = send_idx; i < last_idx; i++)
162 for (i = recv_idx; i < send_idx; i++)
166 Request::sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
168 (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
169 dtype, dst, tag, comm, &status);
171 if (newrank > newdst)
180 // In the non-power-of-two case, all odd-numbered processes of
181 // rank < 2 * rem send the result to (rank-1), the ranks who didn't
182 // participate above.
184 if (rank < 2 * rem) {
186 Request::send(rbuff, count, dtype, rank - 1, tag, comm);
188 Request::recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
191 smpi_free_tmp_buffer(tmp_buf);