Logo AND Algorithmique Numérique Distribuée

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