1 /* Copyright (c) 2013-2014. 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.h"
9 int smpi_coll_tuned_allreduce_rab_rdb(void *sbuff, void *rbuff, int count,
10 MPI_Datatype dtype, MPI_Op op,
13 int tag = COLL_TAG_ALLREDUCE;
14 unsigned int mask, pof2;
15 int dst, newrank, rem, newdst, i,
16 send_idx, recv_idx, last_idx, send_cnt, recv_cnt, *cnts, *disps;
21 unsigned int nprocs = smpi_comm_size(comm);
22 unsigned int rank = smpi_comm_rank(comm);
24 extent = smpi_datatype_get_extent(dtype);
25 tmp_buf = (void *) smpi_get_tmp_sendbuffer(count * extent);
27 smpi_datatype_copy(sbuff, count, dtype, rbuff, count, dtype);
29 // find nearest power-of-two less than or equal to comm_size
31 while (pof2 <= nprocs)
37 // In the non-power-of-two case, all even-numbered
38 // processes of rank < 2*rem send their data to
39 // (rank+1). These even-numbered processes no longer
40 // participate in the algorithm until the very end. The
41 // remaining processes form a nice power-of-two.
47 smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
49 // temporarily set the rank to -1 so that this
50 // process does not pariticipate in recursive
55 smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
56 // do the reduction on received data. since the
57 // ordering is right, it doesn't matter whether
58 // the operation is commutative or not.
59 smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
66 else // rank >= 2 * rem
69 // If op is user-defined or count is less than pof2, use
70 // recursive doubling algorithm. Otherwise do a reduce-scatter
71 // followed by allgather. (If op is user-defined,
72 // derived datatypes are allowed and the user could pass basic
73 // datatypes on one process and derived on another as long as
74 // the type maps are the same. Breaking up derived
75 // datatypes to do the reduce-scatter is tricky, therefore
76 // using recursive doubling in that case.)
79 // do a reduce-scatter followed by allgather. for the
80 // reduce-scatter, calculate the count that each process receives
81 // and the displacement within the buffer
83 cnts = (int *) xbt_malloc(pof2 * sizeof(int));
84 disps = (int *) xbt_malloc(pof2 * sizeof(int));
86 for (i = 0; i < (pof2 - 1); i++)
87 cnts[i] = count / pof2;
88 cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
91 for (i = 1; i < pof2; i++)
92 disps[i] = disps[i - 1] + cnts[i - 1];
95 send_idx = recv_idx = 0;
98 newdst = newrank ^ mask;
99 // find real rank of dest
100 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
102 send_cnt = recv_cnt = 0;
103 if (newrank < newdst) {
104 send_idx = recv_idx + pof2 / (mask * 2);
105 for (i = send_idx; i < last_idx; i++)
107 for (i = recv_idx; i < send_idx; i++)
110 recv_idx = send_idx + pof2 / (mask * 2);
111 for (i = send_idx; i < recv_idx; i++)
113 for (i = recv_idx; i < last_idx; i++)
117 // Send data from recvbuf. Recv into tmp_buf
118 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
120 (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
121 dtype, dst, tag, comm, &status);
123 // tmp_buf contains data received in this step.
124 // recvbuf contains data accumulated so far
126 // This algorithm is used only for predefined ops
127 // and predefined ops are always commutative.
128 smpi_op_apply(op, (char *) tmp_buf + disps[recv_idx] * extent,
129 (char *) rbuff + disps[recv_idx] * extent, &recv_cnt, &dtype);
131 // update send_idx for next iteration
135 // update last_idx, but not in last iteration because the value
136 // is needed in the allgather step below.
138 last_idx = recv_idx + pof2 / mask;
141 // now do the allgather
145 newdst = newrank ^ mask;
146 // find real rank of dest
147 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
149 send_cnt = recv_cnt = 0;
150 if (newrank < newdst) {
151 // update last_idx except on first iteration
152 if (mask != pof2 / 2)
153 last_idx = last_idx + pof2 / (mask * 2);
155 recv_idx = send_idx + pof2 / (mask * 2);
156 for (i = send_idx; i < recv_idx; i++)
158 for (i = recv_idx; i < last_idx; i++)
161 recv_idx = send_idx - pof2 / (mask * 2);
162 for (i = send_idx; i < last_idx; i++)
164 for (i = recv_idx; i < send_idx; i++)
168 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
170 (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
171 dtype, dst, tag, comm, &status);
173 if (newrank > newdst)
183 // In the non-power-of-two case, all odd-numbered processes of
184 // rank < 2 * rem send the result to (rank-1), the ranks who didn't
185 // participate above.
187 if (rank < 2 * rem) {
189 smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
191 smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
194 smpi_free_tmp_buffer(tmp_buf);