Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Remove unused variables.
[simgrid.git] / src / smpi / colls / allreduce-rab-rdb.c
1 #include "colls_private.h"
2
3 int smpi_coll_tuned_allreduce_rab_rdb(void *sbuff, void *rbuff, int count,
4                                       MPI_Datatype dtype, MPI_Op op,
5                                       MPI_Comm comm)
6 {
7   int nprocs, rank, tag = 543;
8   int mask, dst, pof2, newrank, rem, newdst, i,
9       send_idx, recv_idx, last_idx, send_cnt, recv_cnt, *cnts, *disps;
10   MPI_Aint extent;
11   MPI_Status status;
12   void *tmp_buf = NULL;
13
14   nprocs = smpi_comm_size(comm);
15   rank = smpi_comm_rank(comm);
16
17   extent = smpi_datatype_get_extent(dtype);
18   tmp_buf = (void *) xbt_malloc(count * extent);
19
20   smpi_datatype_copy(sbuff, count, dtype, rbuff, count, dtype);
21
22   // find nearest power-of-two less than or equal to comm_size
23   pof2 = 1;
24   while (pof2 <= nprocs)
25     pof2 <<= 1;
26   pof2 >>= 1;
27
28   rem = nprocs - pof2;
29
30   // In the non-power-of-two case, all even-numbered
31   // processes of rank < 2*rem send their data to
32   // (rank+1). These even-numbered processes no longer
33   // participate in the algorithm until the very end. The
34   // remaining processes form a nice power-of-two. 
35
36   if (rank < 2 * rem) {
37     // even       
38     if (rank % 2 == 0) {
39
40       smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
41
42       // temporarily set the rank to -1 so that this
43       // process does not pariticipate in recursive
44       // doubling
45       newrank = -1;
46     } else                      // odd
47     {
48       smpi_mpi_recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
49       // do the reduction on received data. since the
50       // ordering is right, it doesn't matter whether
51       // the operation is commutative or not.
52        smpi_op_apply(op, tmp_buf, rbuff, &count, &dtype);
53
54       // change the rank 
55       newrank = rank / 2;
56     }
57   }
58
59   else                          // rank >= 2 * rem 
60     newrank = rank - rem;
61
62   // If op is user-defined or count is less than pof2, use
63   // recursive doubling algorithm. Otherwise do a reduce-scatter
64   // followed by allgather. (If op is user-defined,
65   // derived datatypes are allowed and the user could pass basic
66   // datatypes on one process and derived on another as long as
67   // the type maps are the same. Breaking up derived
68   // datatypes to do the reduce-scatter is tricky, therefore
69   // using recursive doubling in that case.) 
70
71   if (newrank != -1) {
72     // do a reduce-scatter followed by allgather. for the
73     // reduce-scatter, calculate the count that each process receives
74     // and the displacement within the buffer 
75
76     cnts = (int *) xbt_malloc(pof2 * sizeof(int));
77     disps = (int *) xbt_malloc(pof2 * sizeof(int));
78
79     for (i = 0; i < (pof2 - 1); i++)
80       cnts[i] = count / pof2;
81     cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
82
83     disps[0] = 0;
84     for (i = 1; i < pof2; i++)
85       disps[i] = disps[i - 1] + cnts[i - 1];
86
87     mask = 0x1;
88     send_idx = recv_idx = 0;
89     last_idx = pof2;
90     while (mask < pof2) {
91       newdst = newrank ^ mask;
92       // find real rank of dest 
93       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
94
95       send_cnt = recv_cnt = 0;
96       if (newrank < newdst) {
97         send_idx = recv_idx + pof2 / (mask * 2);
98         for (i = send_idx; i < last_idx; i++)
99           send_cnt += cnts[i];
100         for (i = recv_idx; i < send_idx; i++)
101           recv_cnt += cnts[i];
102       } else {
103         recv_idx = send_idx + pof2 / (mask * 2);
104         for (i = send_idx; i < recv_idx; i++)
105           send_cnt += cnts[i];
106         for (i = recv_idx; i < last_idx; i++)
107           recv_cnt += cnts[i];
108       }
109
110       // Send data from recvbuf. Recv into tmp_buf 
111       smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
112                    dtype, dst, tag,
113                    (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
114                    dtype, dst, tag, comm, &status);
115
116       // tmp_buf contains data received in this step.
117       // recvbuf contains data accumulated so far 
118
119       // This algorithm is used only for predefined ops
120       // and predefined ops are always commutative.
121       smpi_op_apply(op, (char *) tmp_buf + disps[recv_idx] * extent,
122                         (char *) rbuff + disps[recv_idx] * extent, &recv_cnt, &dtype);
123
124       // update send_idx for next iteration 
125       send_idx = recv_idx;
126       mask <<= 1;
127
128       // update last_idx, but not in last iteration because the value
129       // is needed in the allgather step below. 
130       if (mask < pof2)
131         last_idx = recv_idx + pof2 / mask;
132     }
133
134     // now do the allgather 
135
136     mask >>= 1;
137     while (mask > 0) {
138       newdst = newrank ^ mask;
139       // find real rank of dest
140       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
141
142       send_cnt = recv_cnt = 0;
143       if (newrank < newdst) {
144         // update last_idx except on first iteration 
145         if (mask != pof2 / 2)
146           last_idx = last_idx + pof2 / (mask * 2);
147
148         recv_idx = send_idx + pof2 / (mask * 2);
149         for (i = send_idx; i < recv_idx; i++)
150           send_cnt += cnts[i];
151         for (i = recv_idx; i < last_idx; i++)
152           recv_cnt += cnts[i];
153       } else {
154         recv_idx = send_idx - pof2 / (mask * 2);
155         for (i = send_idx; i < last_idx; i++)
156           send_cnt += cnts[i];
157         for (i = recv_idx; i < send_idx; i++)
158           recv_cnt += cnts[i];
159       }
160
161       smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
162                    dtype, dst, tag,
163                    (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
164                    dtype, dst, tag, comm, &status);
165
166       if (newrank > newdst)
167         send_idx = recv_idx;
168
169       mask >>= 1;
170     }
171
172     free(cnts);
173     free(disps);
174
175   }
176   // In the non-power-of-two case, all odd-numbered processes of
177   // rank < 2 * rem send the result to (rank-1), the ranks who didn't
178   // participate above.
179
180   if (rank < 2 * rem) {
181     if (rank % 2)               // odd 
182       smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
183     else                        // even 
184       smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
185   }
186
187   free(tmp_buf);
188   return MPI_SUCCESS;
189 }