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"
15 int Coll_reduce_scatter_gather::reduce(const void *sendbuf, void *recvbuf,
16 int count, MPI_Datatype datatype,
17 MPI_Op op, int root, MPI_Comm comm)
20 int comm_size, rank, pof2, rem, newrank;
21 int mask, *cnts, *disps, i, j, send_idx = 0;
22 int recv_idx, last_idx = 0, newdst;
23 int dst, send_cnt, recv_cnt, newroot, newdst_tree_root;
24 int newroot_tree_root, new_count;
25 int tag = COLL_TAG_REDUCE,temporary_buffer=0;
26 unsigned char *send_ptr, *recv_ptr, *tmp_buf;
36 comm_size = comm->size();
40 extent = datatype->get_extent();
41 /* If I'm not the root, then my recvbuf may not be valid, therefore
42 I have to allocate a temporary one */
43 if (rank != root && not recvbuf) {
45 recvbuf = (void *)smpi_get_tmp_recvbuffer(count * extent);
47 /* find nearest power-of-two less than or equal to comm_size */
49 while (pof2 <= comm_size)
53 if (count < comm_size) {
54 new_count = comm_size;
55 send_ptr = smpi_get_tmp_sendbuffer(new_count * extent);
56 recv_ptr = smpi_get_tmp_recvbuffer(new_count * extent);
57 tmp_buf = smpi_get_tmp_sendbuffer(new_count * extent);
58 memcpy(send_ptr, sendbuf != MPI_IN_PLACE ? sendbuf : recvbuf, extent * count);
61 Request::sendrecv(send_ptr, new_count, datatype, rank, tag,
62 recv_ptr, new_count, datatype, rank, tag, comm, &status);
64 rem = comm_size - pof2;
68 Request::send(recv_ptr, new_count, datatype, rank - 1, tag, comm);
71 Request::recv(tmp_buf, count, datatype, rank + 1, tag, comm, &status);
72 if(op!=MPI_OP_NULL) op->apply( tmp_buf, recv_ptr, &new_count, datatype);
75 } else /* rank >= 2*rem */
79 disps = new int[pof2];
82 for (i = 0; i < (pof2 - 1); i++)
83 cnts[i] = new_count / pof2;
84 cnts[pof2 - 1] = new_count - (new_count / pof2) * (pof2 - 1);
87 for (i = 1; i < pof2; i++)
88 disps[i] = disps[i - 1] + cnts[i - 1];
91 send_idx = recv_idx = 0;
94 newdst = newrank ^ mask;
95 /* find real rank of dest */
96 dst = (newdst < rem) ? newdst * 2 : newdst + rem;
98 send_cnt = recv_cnt = 0;
99 if (newrank < newdst) {
100 send_idx = recv_idx + pof2 / (mask * 2);
101 for (i = send_idx; i < last_idx; i++)
103 for (i = recv_idx; i < send_idx; i++)
106 recv_idx = send_idx + pof2 / (mask * 2);
107 for (i = send_idx; i < recv_idx; i++)
109 for (i = recv_idx; i < last_idx; i++)
113 /* Send data from recvbuf. Recv into tmp_buf */
114 Request::sendrecv(recv_ptr + disps[send_idx] * extent, send_cnt, datatype, dst, tag,
115 tmp_buf + disps[recv_idx] * extent, recv_cnt, datatype, dst, tag, comm, &status);
117 /* tmp_buf contains data received in this step.
118 recvbuf contains data accumulated so far */
120 if (op != MPI_OP_NULL)
121 op->apply(tmp_buf + disps[recv_idx] * extent, recv_ptr + disps[recv_idx] * extent, &recv_cnt, datatype);
123 /* update send_idx for next iteration */
128 last_idx = recv_idx + pof2 / mask;
132 /* now do the gather to root */
134 if (root < 2 * rem) {
138 for (i = 0; i < (pof2 - 1); i++)
139 cnts[i] = new_count / pof2;
140 cnts[pof2 - 1] = new_count - (new_count / pof2) * (pof2 - 1);
143 for (i = 1; i < pof2; i++)
144 disps[i] = disps[i - 1] + cnts[i - 1];
146 Request::recv(recv_ptr, cnts[0], datatype, 0, tag, comm, &status);
151 } else if (newrank == 0) {
152 Request::send(recv_ptr, cnts[0], datatype, root, tag, comm);
159 newroot = root - rem;
164 while (mask < pof2) {
171 newdst = newrank ^ mask;
173 /* find real rank of dest */
174 dst = (newdst < rem) ? newdst * 2 : newdst + rem;
176 if ((newdst == 0) && (root < 2 * rem) && (root % 2 != 0))
178 newdst_tree_root = newdst >> j;
179 newdst_tree_root <<= j;
181 newroot_tree_root = newroot >> j;
182 newroot_tree_root <<= j;
184 send_cnt = recv_cnt = 0;
185 if (newrank < newdst) {
186 /* update last_idx except on first iteration */
187 if (mask != pof2 / 2)
188 last_idx = last_idx + pof2 / (mask * 2);
190 recv_idx = send_idx + pof2 / (mask * 2);
191 for (i = send_idx; i < recv_idx; i++)
193 for (i = recv_idx; i < last_idx; i++)
196 recv_idx = send_idx - pof2 / (mask * 2);
197 for (i = send_idx; i < last_idx; i++)
199 for (i = recv_idx; i < send_idx; i++)
203 if (newdst_tree_root == newroot_tree_root) {
204 Request::send(recv_ptr + disps[send_idx] * extent, send_cnt, datatype, dst, tag, comm);
207 Request::recv(recv_ptr + disps[recv_idx] * extent, recv_cnt, datatype, dst, tag, comm, &status);
210 if (newrank > newdst)
217 memcpy(recvbuf, recv_ptr, extent * count);
218 smpi_free_tmp_buffer(send_ptr);
219 smpi_free_tmp_buffer(recv_ptr);
223 else /* (count >= comm_size) */ {
224 tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
226 //if ((rank != root))
227 Request::sendrecv(sendbuf != MPI_IN_PLACE ? sendbuf : recvbuf, count, datatype, rank, tag,
228 recvbuf, count, datatype, rank, tag, comm, &status);
230 rem = comm_size - pof2;
231 if (rank < 2 * rem) {
232 if (rank % 2 != 0) { /* odd */
233 Request::send(recvbuf, count, datatype, rank - 1, tag, comm);
238 Request::recv(tmp_buf, count, datatype, rank + 1, tag, comm, &status);
239 if(op!=MPI_OP_NULL) op->apply( tmp_buf, recvbuf, &count, datatype);
242 } else /* rank >= 2*rem */
243 newrank = rank - rem;
245 cnts = new int[pof2];
246 disps = new int[pof2];
249 for (i = 0; i < (pof2 - 1); i++)
250 cnts[i] = count / pof2;
251 cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
254 for (i = 1; i < pof2; i++)
255 disps[i] = disps[i - 1] + cnts[i - 1];
258 send_idx = recv_idx = 0;
260 while (mask < pof2) {
261 newdst = newrank ^ mask;
262 /* find real rank of dest */
263 dst = (newdst < rem) ? newdst * 2 : newdst + rem;
265 send_cnt = recv_cnt = 0;
266 if (newrank < newdst) {
267 send_idx = recv_idx + pof2 / (mask * 2);
268 for (i = send_idx; i < last_idx; i++)
270 for (i = recv_idx; i < send_idx; i++)
273 recv_idx = send_idx + pof2 / (mask * 2);
274 for (i = send_idx; i < recv_idx; i++)
276 for (i = recv_idx; i < last_idx; i++)
280 /* Send data from recvbuf. Recv into tmp_buf */
281 Request::sendrecv(static_cast<char*>(recvbuf) + disps[send_idx] * extent, send_cnt, datatype, dst, tag,
282 tmp_buf + disps[recv_idx] * extent, recv_cnt, datatype, dst, tag, comm, &status);
284 /* tmp_buf contains data received in this step.
285 recvbuf contains data accumulated so far */
287 if (op != MPI_OP_NULL)
288 op->apply(tmp_buf + disps[recv_idx] * extent, static_cast<char*>(recvbuf) + disps[recv_idx] * extent,
289 &recv_cnt, datatype);
291 /* update send_idx for next iteration */
296 last_idx = recv_idx + pof2 / mask;
300 /* now do the gather to root */
302 if (root < 2 * rem) {
304 if (rank == root) { /* recv */
305 for (i = 0; i < (pof2 - 1); i++)
306 cnts[i] = count / pof2;
307 cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
310 for (i = 1; i < pof2; i++)
311 disps[i] = disps[i - 1] + cnts[i - 1];
313 Request::recv(recvbuf, cnts[0], datatype, 0, tag, comm, &status);
318 } else if (newrank == 0) {
319 Request::send(recvbuf, cnts[0], datatype, root, tag, comm);
326 newroot = root - rem;
331 while (mask < pof2) {
338 newdst = newrank ^ mask;
340 /* find real rank of dest */
341 dst = (newdst < rem) ? newdst * 2 : newdst + rem;
343 if ((newdst == 0) && (root < 2 * rem) && (root % 2 != 0))
345 newdst_tree_root = newdst >> j;
346 newdst_tree_root <<= j;
348 newroot_tree_root = newroot >> j;
349 newroot_tree_root <<= j;
351 send_cnt = recv_cnt = 0;
352 if (newrank < newdst) {
353 /* update last_idx except on first iteration */
354 if (mask != pof2 / 2)
355 last_idx = last_idx + pof2 / (mask * 2);
357 recv_idx = send_idx + pof2 / (mask * 2);
358 for (i = send_idx; i < recv_idx; i++)
360 for (i = recv_idx; i < last_idx; i++)
363 recv_idx = send_idx - pof2 / (mask * 2);
364 for (i = send_idx; i < last_idx; i++)
366 for (i = recv_idx; i < send_idx; i++)
370 if (newdst_tree_root == newroot_tree_root) {
371 Request::send((char *) recvbuf +
372 disps[send_idx] * extent,
373 send_cnt, datatype, dst, tag, comm);
376 Request::recv((char *) recvbuf +
377 disps[recv_idx] * extent,
378 recv_cnt, datatype, dst, tag, comm, &status);
381 if (newrank > newdst)
390 smpi_free_tmp_buffer(tmp_buf);
391 if (temporary_buffer == 1)
392 smpi_free_tmp_buffer(static_cast<unsigned char*>(recvbuf));