Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge tag 'v3_9_90' into hypervisor
[simgrid.git] / src / smpi / colls / reduce-arrival-pattern-aware.c
1 #include "colls_private.h"
2 //#include <star-reduction.c>
3
4 int reduce_arrival_pattern_aware_segment_size_in_byte = 8192;
5
6 #ifndef HEADER_SIZE
7 #define HEADER_SIZE 1024
8 #endif
9
10 #ifndef MAX_NODE
11 #define MAX_NODE 1024
12 #endif
13
14 /* Non-topology-specific pipelined linear-reduce function */
15 int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf,
16                                                  int count,
17                                                  MPI_Datatype datatype,
18                                                  MPI_Op op, int root,
19                                                  MPI_Comm comm)
20 {
21   int rank = smpi_comm_rank(comm);
22   int tag = -COLL_TAG_REDUCE;
23   MPI_Status status;
24   MPI_Request request;
25   MPI_Request *send_request_array;
26   MPI_Request *recv_request_array;
27   MPI_Status *send_status_array;
28   MPI_Status *recv_status_array;
29
30   MPI_Status temp_status_array[MAX_NODE];
31
32   int size = smpi_comm_size(comm);
33   int i;
34
35   int sent_count;
36   int header_index;
37   int flag_array[MAX_NODE];
38   int already_received[MAX_NODE];
39
40   int header_buf[HEADER_SIZE];
41   char temp_buf[MAX_NODE];
42
43   MPI_Aint extent, lb;
44   smpi_datatype_extent(datatype, &lb, &extent);
45
46   /* source and destination */
47   int to, from;
48
49   /* segment is segment size in number of elements (not bytes) */
50   int segment = reduce_arrival_pattern_aware_segment_size_in_byte / extent;
51
52   /* pipeline length */
53   int pipe_length = count / segment;
54
55   /* use for buffer offset for sending and receiving data = segment size in byte */
56   int increment = segment * extent;
57
58   /* if the input size is not divisible by segment size => 
59      the small remainder will be done with native implementation */
60   int remainder = count % segment;
61
62
63   /* value == 0 means root has not send data (or header) to the node yet */
64   for (i = 0; i < MAX_NODE; i++) {
65     already_received[i] = 0;
66   }
67
68   char *tmp_buf;
69   tmp_buf = (char *) xbt_malloc(count * extent);
70
71   smpi_mpi_sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank,
72                tag, comm, &status);
73
74
75
76   /* when a message is smaller than a block size => no pipeline */
77   if (count <= segment) {
78
79     if (rank == 0) {
80       sent_count = 0;
81
82       while (sent_count < (size - 1)) {
83
84         for (i = 1; i < size; i++) {
85           if (already_received[i] == 0) {
86             smpi_mpi_iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
87                              MPI_STATUSES_IGNORE);
88             simcall_process_sleep(0.0001);
89             }
90         }
91
92         header_index = 0;
93         /* recv 1-byte message */
94         for (i = 0; i < size; i++) {
95           if (i == rank)
96             continue;
97
98           /* 1-byte message arrive */
99           if ((flag_array[i] == 1) && (already_received[i] == 0)) {
100             smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, comm, &status);
101             header_buf[header_index] = i;
102             header_index++;
103             sent_count++;
104
105
106             //printf("root send to %d recv from %d : data = ",to,from);
107             /*
108                for (i=0;i<=header_index;i++) {
109                printf("%d ",header_buf[i]);
110                }
111                printf("\n");
112              */
113             /* will receive in the next step */
114             already_received[i] = 1;
115           }
116         }
117
118         /* send header followed by receive and reduce data */
119         if (header_index != 0) {
120           header_buf[header_index] = -1;
121           to = header_buf[0];
122           from = header_buf[header_index - 1];
123
124           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
125           smpi_mpi_recv(tmp_buf, count, datatype, from, tag, comm, &status);
126           smpi_op_apply(op, tmp_buf, rbuf, &count, &datatype);
127         }
128       }                         /* while loop */
129     }
130
131     /* root */
132     /* non-root */
133     else {
134
135       /* send 1-byte message to root */
136       smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
137
138       /* wait for header and data, forward when required */
139       smpi_mpi_recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
140                &status);
141       //      smpi_mpi_recv(buf,count,datatype,MPI_ANY_SOURCE,tag,comm,&status);
142
143       /* search for where it is */
144       int myordering = 0;
145       while (rank != header_buf[myordering]) {
146         myordering++;
147       }
148
149       /* forward header */
150       if (header_buf[myordering + 1] != -1) {
151           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
152                  tag, comm);
153       }
154       //printf("node %d ordering %d\n",rank,myordering);
155
156       /* receive, reduce, and forward data */
157
158       /* send only */
159       if (myordering == 0) {
160         if (header_buf[myordering + 1] == -1) {
161           to = 0;
162         } else {
163           to = header_buf[myordering + 1];
164         }
165         smpi_mpi_send(rbuf, count, datatype, to, tag, comm);
166       }
167
168       /* recv, reduce, send */
169       else {
170         if (header_buf[myordering + 1] == -1) {
171           to = 0;
172         } else {
173           to = header_buf[myordering + 1];
174         }
175         from = header_buf[myordering - 1];
176         smpi_mpi_recv(tmp_buf, count, datatype, header_buf[myordering - 1], tag,
177                  comm, &status);
178         smpi_op_apply(op, tmp_buf, rbuf, &count, &datatype);
179         smpi_mpi_send(rbuf, count, datatype, to, tag, comm);
180       }
181     }                           /* non-root */
182   }
183   /* pipeline bcast */
184   else {
185     //    printf("node %d start\n",rank);
186
187     send_request_array =
188         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
189     recv_request_array =
190         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
191     send_status_array =
192         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
193     recv_status_array =
194         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
195
196     if (rank == 0) {
197       sent_count = 0;
198
199       int will_send[MAX_NODE];
200       for (i = 0; i < MAX_NODE; i++)
201         will_send[i] = 0;
202
203       /* loop until all data are received (sent) */
204       while (sent_count < (size - 1)) {
205         int k;
206         for (k = 0; k < 1; k++) {
207           for (i = 1; i < size; i++) {
208             //if (i == rank)
209             //continue;
210             if ((already_received[i] == 0) && (will_send[i] == 0)) {
211                 smpi_mpi_iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
212                          &temp_status_array[i]);
213               if (flag_array[i] == 1) {
214                 will_send[i] = 1;
215                 smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
216                          &status);
217                 //printf("recv from %d\n",i);
218                 i = 1;
219               }
220             }
221           }
222         }                       /* end of probing */
223
224         header_index = 0;
225
226         /* recv 1-byte message */
227         for (i = 1; i < size; i++) {
228           //if (i==rank)
229           //continue;
230           /* message arrived in this round (put in the header) */
231           if ((will_send[i] == 1) && (already_received[i] == 0)) {
232             header_buf[header_index] = i;
233             header_index++;
234             sent_count++;
235
236             /* will send in the next step */
237             already_received[i] = 1;
238           }
239         }
240
241         /* send header followed by data */
242         if (header_index != 0) {
243           header_buf[header_index] = -1;
244           to = header_buf[0];
245
246           /* send header */
247           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
248
249           /* recv data - pipeline */
250           from = header_buf[header_index - 1];
251           for (i = 0; i < pipe_length; i++) {
252             smpi_mpi_recv(tmp_buf + (i * increment), segment, datatype, from, tag,
253                      comm, &status);
254             smpi_op_apply(op, tmp_buf + (i * increment),
255                            (char *)rbuf + (i * increment), &segment, &datatype);
256           }
257         }
258       }                         /* while loop (sent_count < size-1 ) */
259     }
260
261     /* root */
262     /* none root */
263     else {
264       /* send 1-byte message to root */
265       smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
266
267
268       /* wait for header forward when required */
269       request=smpi_mpi_irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm);
270       smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
271
272       /* search for where it is */
273       int myordering = 0;
274
275       while (rank != header_buf[myordering]) {
276         myordering++;
277       }
278
279       /* send header when required */
280       if (header_buf[myordering + 1] != -1) {
281           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
282                  tag, comm);
283       }
284
285       /* (receive, reduce), and send data */
286       if (header_buf[myordering + 1] == -1) {
287         to = 0;
288       } else {
289         to = header_buf[myordering + 1];
290       }
291
292       /* send only */
293       if (myordering == 0) {
294         for (i = 0; i < pipe_length; i++) {
295             send_request_array[i]= smpi_mpi_isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm);
296         }
297         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
298       }
299
300       /* receive, reduce, and send */
301       else {
302         from = header_buf[myordering - 1];
303         for (i = 0; i < pipe_length; i++) {
304           recv_request_array[i]=smpi_mpi_irecv(tmp_buf + (i * increment), segment, datatype, from, tag, comm);
305         }
306         for (i = 0; i < pipe_length; i++) {
307           smpi_mpi_wait(&recv_request_array[i], MPI_STATUS_IGNORE);
308           smpi_op_apply(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment),
309                          &segment, &datatype);
310           send_request_array[i]=smpi_mpi_isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm);
311         }
312         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
313       }
314     }                           /* non-root */
315
316
317
318
319     free(send_request_array);
320     free(recv_request_array);
321     free(send_status_array);
322     free(recv_status_array);
323
324     //printf("node %d done\n",rank);
325   }                             /* end pipeline */
326
327
328   /* if root is not zero send root after finished
329      this can be modified to make it faster by using logical src, dst.
330    */
331   if (root != 0) {
332     if (rank == 0) {
333       smpi_mpi_send(rbuf, count, datatype, root, tag, comm);
334     } else if (rank == root) {
335       smpi_mpi_recv(rbuf, count, datatype, 0, tag, comm, &status);
336     }
337   }
338
339
340   /* when count is not divisible by block size, use default BCAST for the remainder */
341   if ((remainder != 0) && (count > segment)) {
342     smpi_mpi_reduce((char *)buf + (pipe_length * increment),
343                (char *)rbuf + (pipe_length * increment), remainder, datatype, op, root,
344                comm);
345   }
346
347   free(tmp_buf);
348
349   return MPI_SUCCESS;
350 }