1 #include "colls_private.h"
2 //#include <star-reduction.c>
4 int smpi_coll_tuned_allreduce_rab_rsag(void *sbuff, void *rbuff, int count,
5 MPI_Datatype dtype, MPI_Op op,
8 int nprocs, rank, tag = 543;
9 int mask, dst, pof2, newrank, rem, newdst, i,
10 send_idx, recv_idx, last_idx, send_cnt, recv_cnt, *cnts, *disps;
14 nprocs = smpi_comm_size(comm);
15 rank = smpi_comm_rank(comm);
17 extent = smpi_datatype_get_extent(dtype);
18 tmp_buf = (void *) xbt_malloc(count * extent);
20 smpi_mpi_sendrecv(sbuff, count, dtype, rank, tag, rbuff, count, dtype, rank, tag,
23 // find nearest power-of-two less than or equal to comm_size
25 while (pof2 <= nprocs)
31 // In the non-power-of-two case, all even-numbered
32 // processes of rank < 2*rem send their data to
33 // (rank+1). These even-numbered processes no longer
34 // participate in the algorithm until the very end. The
35 // remaining processes form a nice power-of-two.
41 smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
43 // temporarily set the rank to -1 so that this
44 // process does not pariticipate in recursive
49 smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
50 // do the reduction on received data. since the
51 // ordering is right, it doesn't matter whether
52 // the operation is commutative or not.
53 smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
60 else // rank >= 2 * rem
63 // If op is user-defined or count is less than pof2, use
64 // recursive doubling algorithm. Otherwise do a reduce-scatter
65 // followed by allgather. (If op is user-defined,
66 // derived datatypes are allowed and the user could pass basic
67 // datatypes on one process and derived on another as long as
68 // the type maps are the same. Breaking up derived
69 // datatypes to do the reduce-scatter is tricky, therefore
70 // using recursive doubling in that case.)
73 // do a reduce-scatter followed by allgather. for the
74 // reduce-scatter, calculate the count that each process receives
75 // and the displacement within the buffer
77 cnts = (int *) xbt_malloc(pof2 * sizeof(int));
78 disps = (int *) xbt_malloc(pof2 * sizeof(int));
80 for (i = 0; i < (pof2 - 1); i++)
81 cnts[i] = count / pof2;
82 cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
85 for (i = 1; i < pof2; i++)
86 disps[i] = disps[i - 1] + cnts[i - 1];
89 send_idx = recv_idx = 0;
92 newdst = newrank ^ mask;
93 // find real rank of dest
94 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
96 send_cnt = recv_cnt = 0;
97 if (newrank < newdst) {
98 send_idx = recv_idx + pof2 / (mask * 2);
99 for (i = send_idx; i < last_idx; i++)
101 for (i = recv_idx; i < send_idx; i++)
104 recv_idx = send_idx + pof2 / (mask * 2);
105 for (i = send_idx; i < recv_idx; i++)
107 for (i = recv_idx; i < last_idx; i++)
111 // Send data from recvbuf. Recv into tmp_buf
112 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
114 (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
115 dtype, dst, tag, comm, &status);
117 // tmp_buf contains data received in this step.
118 // recvbuf contains data accumulated so far
120 // This algorithm is used only for predefined ops
121 // and predefined ops are always commutative.
122 smpi_op_apply(op, (char *) tmp_buf + disps[recv_idx] * extent,
123 (char *) rbuff + disps[recv_idx] * extent,
126 // update send_idx for next iteration
130 // update last_idx, but not in last iteration because the value
131 // is needed in the allgather step below.
133 last_idx = recv_idx + pof2 / mask;
136 // now do the allgather
140 newdst = newrank ^ mask;
141 // find real rank of dest
142 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
144 send_cnt = recv_cnt = 0;
145 if (newrank < newdst) {
146 // update last_idx except on first iteration
147 if (mask != pof2 / 2)
148 last_idx = last_idx + pof2 / (mask * 2);
150 recv_idx = send_idx + pof2 / (mask * 2);
151 for (i = send_idx; i < recv_idx; i++)
153 for (i = recv_idx; i < last_idx; i++)
156 recv_idx = send_idx - pof2 / (mask * 2);
157 for (i = send_idx; i < last_idx; i++)
159 for (i = recv_idx; i < send_idx; i++)
163 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
165 (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
166 dtype, dst, tag, comm, &status);
168 if (newrank > newdst)
178 // In the non-power-of-two case, all odd-numbered processes of
179 // rank < 2 * rem send the result to (rank-1), the ranks who didn't
180 // participate above.
182 if (rank < 2 * rem) {
184 smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
186 smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);