Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
54149d4db4dfc6bfba58637652cd25c9aa57bc93
[simgrid.git] / src / smpi / colls / allreduce-rab-rsag.c
1 #include "colls.h"
2 //#include <star-reduction.c>
3
4 int smpi_coll_tuned_allreduce_rab_rsag(void *sbuff, void *rbuff, int count,
5                                        MPI_Datatype dtype, MPI_Op op,
6                                        MPI_Comm comm)
7 {
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;
11   MPI_Aint extent;
12   MPI_Status status;
13   void *tmp_buf = NULL;
14   MPI_Comm_size(comm, &nprocs);
15   MPI_Comm_rank(comm, &rank);
16
17   MPI_Type_extent(dtype, &extent);
18   tmp_buf = (void *) xbt_malloc(count * extent);
19
20   MPI_Sendrecv(sbuff, count, dtype, rank, tag, rbuff, count, dtype, rank, tag,
21                comm, &status);
22
23   MPI_Type_size(dtype, &type_size);
24
25   // find nearest power-of-two less than or equal to comm_size
26   pof2 = 1;
27   while (pof2 <= nprocs)
28     pof2 <<= 1;
29   pof2 >>= 1;
30
31   rem = nprocs - pof2;
32
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. 
38
39   if (rank < 2 * rem) {
40     // even       
41     if (rank % 2 == 0) {
42
43       MPI_Send(rbuff, count, dtype, rank + 1, tag, comm);
44
45       // temporarily set the rank to -1 so that this
46       // process does not pariticipate in recursive
47       // doubling
48       newrank = -1;
49     } else                      // odd
50     {
51       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       star_reduction(op, tmp_buf, rbuff, &count, &dtype);
56
57       // change the rank 
58       newrank = rank / 2;
59     }
60   }
61
62   else                          // rank >= 2 * rem 
63     newrank = rank - rem;
64
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.) 
73
74   if (newrank != -1) {
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 
78
79     cnts = (int *) malloc(pof2 * sizeof(int));
80     disps = (int *) malloc(pof2 * sizeof(int));
81
82     for (i = 0; i < (pof2 - 1); i++)
83       cnts[i] = count / pof2;
84     cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
85
86     disps[0] = 0;
87     for (i = 1; i < pof2; i++)
88       disps[i] = disps[i - 1] + cnts[i - 1];
89
90     mask = 0x1;
91     send_idx = recv_idx = 0;
92     last_idx = pof2;
93     while (mask < pof2) {
94       newdst = newrank ^ mask;
95       // find real rank of dest 
96       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
97
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++)
102           send_cnt += cnts[i];
103         for (i = recv_idx; i < send_idx; i++)
104           recv_cnt += cnts[i];
105       } else {
106         recv_idx = send_idx + pof2 / (mask * 2);
107         for (i = send_idx; i < recv_idx; i++)
108           send_cnt += cnts[i];
109         for (i = recv_idx; i < last_idx; i++)
110           recv_cnt += cnts[i];
111       }
112
113       // Send data from recvbuf. Recv into tmp_buf 
114       MPI_Sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
115                    dtype, dst, tag,
116                    (char *) tmp_buf + disps[recv_idx] * extent, recv_cnt,
117                    dtype, dst, tag, comm, &status);
118
119       // tmp_buf contains data received in this step.
120       // recvbuf contains data accumulated so far 
121
122       // This algorithm is used only for predefined ops
123       // and predefined ops are always commutative.
124       star_reduction(op, (char *) tmp_buf + disps[recv_idx] * extent,
125                      (char *) rbuff + disps[recv_idx] * extent,
126                      &recv_cnt, &dtype);
127
128       // update send_idx for next iteration 
129       send_idx = recv_idx;
130       mask <<= 1;
131
132       // update last_idx, but not in last iteration because the value
133       // is needed in the allgather step below. 
134       if (mask < pof2)
135         last_idx = recv_idx + pof2 / mask;
136     }
137
138     // now do the allgather 
139
140     mask >>= 1;
141     while (mask > 0) {
142       newdst = newrank ^ mask;
143       // find real rank of dest
144       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
145
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);
151
152         recv_idx = send_idx + pof2 / (mask * 2);
153         for (i = send_idx; i < recv_idx; i++)
154           send_cnt += cnts[i];
155         for (i = recv_idx; i < last_idx; i++)
156           recv_cnt += cnts[i];
157       } else {
158         recv_idx = send_idx - pof2 / (mask * 2);
159         for (i = send_idx; i < last_idx; i++)
160           send_cnt += cnts[i];
161         for (i = recv_idx; i < send_idx; i++)
162           recv_cnt += cnts[i];
163       }
164
165       MPI_Sendrecv((char *) rbuff + disps[send_idx] * extent, send_cnt,
166                    dtype, dst, tag,
167                    (char *) rbuff + disps[recv_idx] * extent, recv_cnt,
168                    dtype, dst, tag, comm, &status);
169
170       if (newrank > newdst)
171         send_idx = recv_idx;
172
173       mask >>= 1;
174     }
175
176     free(cnts);
177     free(disps);
178
179   }
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.
183
184   if (rank < 2 * rem) {
185     if (rank % 2)               // odd 
186       MPI_Send(rbuff, count, dtype, rank - 1, tag, comm);
187     else                        // even 
188       MPI_Recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
189   }
190
191   free(tmp_buf);
192   return MPI_SUCCESS;
193 }