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 nprocs, rank, tag = COLL_TAG_ALLREDUCE;
14 int mask, dst, pof2, newrank, rem, newdst, i,
15 send_idx, recv_idx, last_idx, send_cnt, recv_cnt, *cnts, *disps;
20 nprocs = smpi_comm_size(comm);
21 rank = smpi_comm_rank(comm);
23 extent = smpi_datatype_get_extent(dtype);
24 tmp_buf = (void *) xbt_malloc(count * extent);
26 smpi_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 smpi_mpi_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 smpi_mpi_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 smpi_op_apply(op, 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 cnts = (int *) xbt_malloc(pof2 * sizeof(int));
83 disps = (int *) xbt_malloc(pof2 * sizeof(int));
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 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
119 (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
120 dtype, dst, tag, comm, &status);
122 // tmp_buf contains data received in this step.
123 // recvbuf contains data accumulated so far
125 // This algorithm is used only for predefined ops
126 // and predefined ops are always commutative.
127 smpi_op_apply(op, (char *) tmp_buf + disps[recv_idx] * extent,
128 (char *) rbuff + disps[recv_idx] * extent, &recv_cnt, &dtype);
130 // update send_idx for next iteration
134 // update last_idx, but not in last iteration because the value
135 // is needed in the allgather step below.
137 last_idx = recv_idx + pof2 / mask;
140 // now do the allgather
144 newdst = newrank ^ mask;
145 // find real rank of dest
146 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
148 send_cnt = recv_cnt = 0;
149 if (newrank < newdst) {
150 // update last_idx except on first iteration
151 if (mask != pof2 / 2)
152 last_idx = last_idx + pof2 / (mask * 2);
154 recv_idx = send_idx + pof2 / (mask * 2);
155 for (i = send_idx; i < recv_idx; i++)
157 for (i = recv_idx; i < last_idx; i++)
160 recv_idx = send_idx - pof2 / (mask * 2);
161 for (i = send_idx; i < last_idx; i++)
163 for (i = recv_idx; i < send_idx; i++)
167 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
169 (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
170 dtype, dst, tag, comm, &status);
172 if (newrank > newdst)
182 // In the non-power-of-two case, all odd-numbered processes of
183 // rank < 2 * rem send the result to (rank-1), the ranks who didn't
184 // participate above.
186 if (rank < 2 * rem) {
188 smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
190 smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);