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, type_size, 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 type_size = smpi_datatype_size(dtype);
25 // find nearest power-of-two less than or equal to comm_size
27 while (pof2 <= nprocs)
33 // In the non-power-of-two case, all even-numbered
34 // processes of rank < 2*rem send their data to
35 // (rank+1). These even-numbered processes no longer
36 // participate in the algorithm until the very end. The
37 // remaining processes form a nice power-of-two.
43 smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
45 // temporarily set the rank to -1 so that this
46 // process does not pariticipate in recursive
51 smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
52 // do the reduction on received data. since the
53 // ordering is right, it doesn't matter whether
54 // the operation is commutative or not.
55 smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
62 else // rank >= 2 * rem
65 // If op is user-defined or count is less than pof2, use
66 // recursive doubling algorithm. Otherwise do a reduce-scatter
67 // followed by allgather. (If op is user-defined,
68 // derived datatypes are allowed and the user could pass basic
69 // datatypes on one process and derived on another as long as
70 // the type maps are the same. Breaking up derived
71 // datatypes to do the reduce-scatter is tricky, therefore
72 // using recursive doubling in that case.)
75 // do a reduce-scatter followed by allgather. for the
76 // reduce-scatter, calculate the count that each process receives
77 // and the displacement within the buffer
79 cnts = (int *) xbt_malloc(pof2 * sizeof(int));
80 disps = (int *) xbt_malloc(pof2 * sizeof(int));
82 for (i = 0; i < (pof2 - 1); i++)
83 cnts[i] = count / pof2;
84 cnts[pof2 - 1] = count - (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 + 1 : 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 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
116 (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
117 dtype, dst, tag, comm, &status);
119 // tmp_buf contains data received in this step.
120 // recvbuf contains data accumulated so far
122 // This algorithm is used only for predefined ops
123 // and predefined ops are always commutative.
124 smpi_op_apply(op, (char *) tmp_buf + disps[recv_idx] * extent,
125 (char *) rbuff + disps[recv_idx] * extent,
128 // update send_idx for next iteration
132 // update last_idx, but not in last iteration because the value
133 // is needed in the allgather step below.
135 last_idx = recv_idx + pof2 / mask;
138 // now do the allgather
142 newdst = newrank ^ mask;
143 // find real rank of dest
144 dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
146 send_cnt = recv_cnt = 0;
147 if (newrank < newdst) {
148 // update last_idx except on first iteration
149 if (mask != pof2 / 2)
150 last_idx = last_idx + pof2 / (mask * 2);
152 recv_idx = send_idx + pof2 / (mask * 2);
153 for (i = send_idx; i < recv_idx; i++)
155 for (i = recv_idx; i < last_idx; i++)
158 recv_idx = send_idx - pof2 / (mask * 2);
159 for (i = send_idx; i < last_idx; i++)
161 for (i = recv_idx; i < send_idx; i++)
165 smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
167 (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
168 dtype, dst, tag, comm, &status);
170 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 smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
188 smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);