Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
3b0a60595ad5dcb5a01a65c7bf5bce492fe458ef
[simgrid.git] / src / smpi / colls / allreduce / allreduce-rab-rdb.cpp
1 /* Copyright (c) 2013-2022. 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.hpp"
8 namespace simgrid::smpi {
9 int allreduce__rab_rdb(const 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, i, recv_idx, last_idx, send_idx, send_cnt;
15   int dst, newrank, rem, newdst, recv_cnt;
16   MPI_Aint extent;
17   MPI_Status status;
18
19   unsigned int nprocs = comm->size();
20   int rank = comm->rank();
21
22   extent = dtype->get_extent();
23   unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
24
25   Datatype::copy(sbuff, count, dtype, rbuff, count, dtype);
26
27   // find nearest power-of-two less than or equal to comm_size
28   pof2 = 1;
29   while (pof2 <= nprocs)
30     pof2 <<= 1;
31   pof2 >>= 1;
32
33   rem = nprocs - pof2;
34
35   // In the non-power-of-two case, all even-numbered
36   // processes of rank < 2*rem send their data to
37   // (rank+1). These even-numbered processes no longer
38   // participate in the algorithm until the very end. The
39   // remaining processes form a nice power-of-two.
40
41   if (rank < 2 * rem) {
42     // even
43     if (rank % 2 == 0) {
44
45       Request::send(rbuff, count, dtype, rank + 1, tag, comm);
46
47       // temporarily set the rank to -1 so that this
48       // process does not participate in recursive
49       // doubling
50       newrank = -1;
51     } else                      // odd
52     {
53       Request::recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
54       // do the reduction on received data. since the
55       // ordering is right, it doesn't matter whether
56       // the operation is commutative or not.
57        if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuff, &count, dtype);
58
59       // change the rank
60       newrank = rank / 2;
61     }
62   }
63
64   else                          // rank >= 2 * rem
65     newrank = rank - rem;
66
67   // If op is user-defined or count is less than pof2, use
68   // recursive doubling algorithm. Otherwise do a reduce-scatter
69   // followed by allgather. (If op is user-defined,
70   // derived datatypes are allowed and the user could pass basic
71   // datatypes on one process and derived on another as long as
72   // the type maps are the same. Breaking up derived
73   // datatypes to do the reduce-scatter is tricky, therefore
74   // using recursive doubling in that case.)
75
76   if (newrank != -1) {
77     // do a reduce-scatter followed by allgather. for the
78     // reduce-scatter, calculate the count that each process receives
79     // and the displacement within the buffer
80
81     int* cnts  = new int[pof2];
82     int* disps = new int[pof2];
83
84     for (i = 0; i < (pof2 - 1); i++)
85       cnts[i] = count / pof2;
86     cnts[pof2 - 1] = count - (count / pof2) * (pof2 - 1);
87
88     disps[0] = 0;
89     for (i = 1; i < pof2; i++)
90       disps[i] = disps[i - 1] + cnts[i - 1];
91
92     mask = 0x1;
93     send_idx = recv_idx = 0;
94     last_idx = pof2;
95     while (mask < pof2) {
96       newdst = newrank ^ mask;
97       // find real rank of dest
98       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
99
100       send_cnt = recv_cnt = 0;
101       if (newrank < newdst) {
102         send_idx = recv_idx + pof2 / (mask * 2);
103         for (i = send_idx; i < last_idx; i++)
104           send_cnt += cnts[i];
105         for (i = recv_idx; i < send_idx; i++)
106           recv_cnt += cnts[i];
107       } else {
108         recv_idx = send_idx + pof2 / (mask * 2);
109         for (i = send_idx; i < recv_idx; i++)
110           send_cnt += cnts[i];
111         for (i = recv_idx; i < last_idx; i++)
112           recv_cnt += cnts[i];
113       }
114
115       // Send data from recvbuf. Recv into tmp_buf
116       Request::sendrecv(static_cast<char*>(rbuff) + disps[send_idx] * extent, send_cnt, dtype, dst, tag,
117                         tmp_buf + disps[recv_idx] * extent, recv_cnt, 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       if (op != MPI_OP_NULL)
125         op->apply(tmp_buf + disps[recv_idx] * extent, static_cast<char*>(rbuff) + disps[recv_idx] * extent, &recv_cnt,
126                   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       Request::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     delete[] cnts;
177     delete[] disps;
178   }
179   // In the non-power-of-two case, all odd-numbered processes of
180   // rank < 2 * rem send the result to (rank-1), the ranks who didn't
181   // participate above.
182
183   if (rank < 2 * rem) {
184     if (rank % 2)               // odd
185       Request::send(rbuff, count, dtype, rank - 1, tag, comm);
186     else                        // even
187       Request::recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
188   }
189
190   smpi_free_tmp_buffer(tmp_buf);
191   return MPI_SUCCESS;
192 }
193 } // namespace simgrid::smpi