Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
7499785afb09205d9f821f087a849a14ba5d7299
[simgrid.git] / src / smpi / colls / allreduce-rab-rdb.c
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
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. */
6
7 #include "colls_private.h"
8
9 int smpi_coll_tuned_allreduce_rab_rdb(void *sbuff, void *rbuff, int count,
10                                       MPI_Datatype dtype, MPI_Op op,
11                                       MPI_Comm comm)
12 {
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;
17   MPI_Aint extent;
18   MPI_Status status;
19   void *tmp_buf = NULL;
20
21   unsigned int nprocs = smpi_comm_size(comm);
22   unsigned int rank = smpi_comm_rank(comm);
23
24   extent = smpi_datatype_get_extent(dtype);
25   tmp_buf = (void *) smpi_get_tmp_sendbuffer(count * extent);
26
27   smpi_datatype_copy(sbuff, count, dtype, rbuff, count, dtype);
28
29   // find nearest power-of-two less than or equal to comm_size
30   pof2 = 1;
31   while (pof2 <= nprocs)
32     pof2 <<= 1;
33   pof2 >>= 1;
34
35   rem = nprocs - pof2;
36
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. 
42
43   if (rank < 2 * rem) {
44     // even       
45     if (rank % 2 == 0) {
46
47       smpi_mpi_send(rbuff, count, dtype, rank + 1, tag, comm);
48
49       // temporarily set the rank to -1 so that this
50       // process does not pariticipate in recursive
51       // doubling
52       newrank = -1;
53     } else                      // odd
54     {
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);
60
61       // change the rank 
62       newrank = rank / 2;
63     }
64   }
65
66   else                          // rank >= 2 * rem 
67     newrank = rank - rem;
68
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.) 
77
78   if (newrank != -1) {
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 
82
83     cnts = (int *) xbt_malloc(pof2 * sizeof(int));
84     disps = (int *) xbt_malloc(pof2 * sizeof(int));
85
86     for (i = 0; i < (pof2 - 1); i++)
87       cnts[i] = count / pof2;
88     cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
89
90     disps[0] = 0;
91     for (i = 1; i < pof2; i++)
92       disps[i] = disps[i - 1] + cnts[i - 1];
93
94     mask = 0x1;
95     send_idx = recv_idx = 0;
96     last_idx = pof2;
97     while (mask < pof2) {
98       newdst = newrank ^ mask;
99       // find real rank of dest 
100       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
101
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++)
106           send_cnt += cnts[i];
107         for (i = recv_idx; i < send_idx; i++)
108           recv_cnt += cnts[i];
109       } else {
110         recv_idx = send_idx + pof2 / (mask * 2);
111         for (i = send_idx; i < recv_idx; i++)
112           send_cnt += cnts[i];
113         for (i = recv_idx; i < last_idx; i++)
114           recv_cnt += cnts[i];
115       }
116
117       // Send data from recvbuf. Recv into tmp_buf 
118       smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
119                    dtype, dst, tag,
120                    (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
121                    dtype, dst, tag, comm, &status);
122
123       // tmp_buf contains data received in this step.
124       // recvbuf contains data accumulated so far 
125
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);
130
131       // update send_idx for next iteration 
132       send_idx = recv_idx;
133       mask <<= 1;
134
135       // update last_idx, but not in last iteration because the value
136       // is needed in the allgather step below. 
137       if (mask < pof2)
138         last_idx = recv_idx + pof2 / mask;
139     }
140
141     // now do the allgather 
142
143     mask >>= 1;
144     while (mask > 0) {
145       newdst = newrank ^ mask;
146       // find real rank of dest
147       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
148
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);
154
155         recv_idx = send_idx + pof2 / (mask * 2);
156         for (i = send_idx; i < recv_idx; i++)
157           send_cnt += cnts[i];
158         for (i = recv_idx; i < last_idx; i++)
159           recv_cnt += cnts[i];
160       } else {
161         recv_idx = send_idx - pof2 / (mask * 2);
162         for (i = send_idx; i < last_idx; i++)
163           send_cnt += cnts[i];
164         for (i = recv_idx; i < send_idx; i++)
165           recv_cnt += cnts[i];
166       }
167
168       smpi_mpi_sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
169                    dtype, dst, tag,
170                    (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
171                    dtype, dst, tag, comm, &status);
172
173       if (newrank > newdst)
174         send_idx = recv_idx;
175
176       mask >>= 1;
177     }
178
179     free(cnts);
180     free(disps);
181
182   }
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.
186
187   if (rank < 2 * rem) {
188     if (rank % 2)               // odd 
189       smpi_mpi_send(rbuff, count, dtype, rank - 1, tag, comm);
190     else                        // even 
191       smpi_mpi_recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
192   }
193
194   smpi_free_tmp_buffer(tmp_buf);
195   return MPI_SUCCESS;
196 }