Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Make internal collective flags negative (incorrect in MPI), to avoid confusion with...
[simgrid.git] / src / smpi / colls / bcast-arrival-nb.c
1 #include "colls_private.h"
2
3 static int bcast_NTSL_segment_size_in_byte = 8192;
4
5 #define HEADER_SIZE 1024
6 #define MAX_NODE 1024
7
8 /* Non-topology-specific pipelined linear-bcast function */
9 int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
10                                      MPI_Datatype datatype, int root,
11                                      MPI_Comm comm)
12 {
13   int tag = -COLL_TAG_BCAST;
14   MPI_Status status;
15   MPI_Request request;
16   MPI_Request *send_request_array;
17   MPI_Request *recv_request_array;
18   MPI_Status *send_status_array;
19   MPI_Status *recv_status_array;
20
21   MPI_Status temp_status_array[MAX_NODE];
22
23   int rank, size;
24   int i, j;
25
26   int sent_count;
27   int header_index;
28   int flag_array[MAX_NODE];
29   int already_sent[MAX_NODE];
30   int to_clean[MAX_NODE];
31   int header_buf[HEADER_SIZE];
32   char temp_buf[MAX_NODE];
33
34   MPI_Aint extent;
35   extent = smpi_datatype_get_extent(datatype);
36
37   /* destination */
38   int to;
39
40
41
42   rank = smpi_comm_rank(MPI_COMM_WORLD);
43   size = smpi_comm_size(MPI_COMM_WORLD);
44
45
46   /* segment is segment size in number of elements (not bytes) */
47   int segment = bcast_NTSL_segment_size_in_byte / extent;
48
49   /* pipeline length */
50   int pipe_length = count / segment;
51
52   /* use for buffer offset for sending and receiving data = segment size in byte */
53   int increment = segment * extent;
54
55   /* if the input size is not divisible by segment size => 
56      the small remainder will be done with native implementation */
57   int remainder = count % segment;
58
59   /* if root is not zero send to rank zero first
60      this can be modified to make it faster by using logical src, dst.
61    */
62   if (root != 0) {
63     if (rank == root) {
64       smpi_mpi_send(buf, count, datatype, 0, tag, comm);
65     } else if (rank == 0) {
66       smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
67     }
68   }
69
70   /* value == 0 means root has not send data (or header) to the node yet */
71   for (i = 0; i < MAX_NODE; i++) {
72     already_sent[i] = 0;
73     to_clean[i]=0;
74   }
75   //  printf("YYY\n");
76
77   /* when a message is smaller than a block size => no pipeline */
78   if (count <= segment) {
79     if (rank == 0) {
80       sent_count = 0;
81
82       while (sent_count < (size - 1)) {
83
84         //      for (j=0;j<1000;j++) {
85         for (i = 1; i < size; i++) {
86           if (already_sent[i] == 0)
87             smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
88                        MPI_STATUSES_IGNORE);
89         }
90         //}
91
92         header_index = 0;
93         /* recv 1-byte message */
94         for (i = 1; i < size; i++) {
95
96           /* message arrive */
97           if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
98             smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
99             header_buf[header_index] = i;
100             header_index++;
101             sent_count++;
102
103             /* will send in the next step */
104             already_sent[i] = 1;
105           }
106         }
107
108         /* send header followed by data */
109         if (header_index != 0) {
110           header_buf[header_index] = -1;
111           to = header_buf[0];
112           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
113           smpi_mpi_send(buf, count, datatype, to, tag, comm);
114         }
115
116         /* randomly MPI_Send to one */
117         else {
118           /* search for the first node that never received data before */
119           for (i = 1; i < size; i++) {
120             if (already_sent[i] == 0) {
121               header_buf[0] = i;
122               header_buf[1] = -1;
123               smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
124               smpi_mpi_send(buf, count, datatype, i, tag, comm);
125               already_sent[i] = 1;
126               sent_count++;
127               to_clean[i]=0;
128               break;
129             }
130           }
131         }
132
133
134       }                         /* while loop */
135       
136       for(i=0; i<size; i++)
137         if(to_clean[i]!=0)smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
138     }
139
140     /* non-root */
141     else {
142
143       /* send 1-byte message to root */
144       smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
145
146       /* wait for header and data, forward when required */
147       smpi_mpi_recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
148                &status);
149       smpi_mpi_recv(buf, count, datatype, MPI_ANY_SOURCE, tag, comm, &status);
150
151       /* search for where it is */
152       int myordering = 0;
153       while (rank != header_buf[myordering]) {
154         myordering++;
155       }
156
157       /* send header followed by data */
158       if (header_buf[myordering + 1] != -1) {
159         smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
160                  tag, comm);
161         smpi_mpi_send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
162       }
163     }
164   }
165   /* pipeline bcast */
166   else {
167     send_request_array =
168         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
169     recv_request_array =
170         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
171     send_status_array =
172         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
173     recv_status_array =
174         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
175
176     if (rank == 0) {
177       sent_count = 0;
178       int iteration = 0;
179
180       int will_send[1000];
181       for (i = 0; i < 1000; i++)
182         will_send[i] = 0;
183       while (sent_count < (size - 1)) {
184         iteration++;
185         //start = MPI_Wtime();
186
187         int k;
188         for (k = 0; k < 3; k++) {
189           for (i = 1; i < size; i++) {
190             if ((already_sent[i] == 0) && (will_send[i] == 0)) {
191               smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
192                          &temp_status_array[i]);
193               if (flag_array[i] == 1) {
194                 will_send[i] = 1;
195                 smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD,
196                          &status);
197                 i = 1;
198               }
199             }
200           }
201         }
202
203         //total = MPI_Wtime() - start;
204         //total *= 1000;
205         //printf("Iprobe time = %.2f\n",total);
206         header_index = 0;
207
208         //start = MPI_Wtime();
209         /* recv 1-byte message */
210         for (i = 1; i < size; i++) {
211           /* message arrive */
212           if ((will_send[i] == 1) && (already_sent[i] == 0)) {
213             header_buf[header_index] = i;
214             header_index++;
215             sent_count++;
216
217             /* will send in the next step */
218             already_sent[i] = 1;
219           }
220         }
221         //printf("sent_count = %d\n",sent_count);
222
223
224         //total = MPI_Wtime() - start;
225         //total *= 1000;
226         //printf("Recv 1-byte time = %.2f\n",total);
227
228         /*      
229            if (header_index != 0) {
230            printf("header index = %d node = ",header_index);
231            for (i=0;i<header_index;i++) {
232            printf("%d ",header_buf[i]);
233            }
234            printf("\n");
235            }
236          */
237
238         /* send header followed by data */
239         if (header_index != 0) {
240           header_buf[header_index] = -1;
241           to = header_buf[0];
242
243           //start = MPI_Wtime();
244
245           /* send header */
246           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
247
248           //total = MPI_Wtime() - start;
249           //total *= 1000;
250           //printf("\tSend header to %d time = %.2f\n",to,total);
251
252           //start = MPI_Wtime();
253
254           /* send data - non-pipeline case */
255
256           if (0 == 1) {
257             //if (header_index == 1) {
258             smpi_mpi_send(buf, count, datatype, to, tag, comm);
259           }
260
261
262           /* send data - pipeline */
263           else {
264             for (i = 0; i < pipe_length; i++) {
265               smpi_mpi_send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
266             }
267             //smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
268           }
269           //total = MPI_Wtime() - start;
270           //total *= 1000;
271           //printf("\tSend data to %d time = %.2f\n",to,total);
272
273         }
274
275
276
277         /* randomly MPI_Send to one node */
278         else {
279           /* search for the first node that never received data before */
280           for (i = 1; i < size; i++) {
281             if (already_sent[i] == 0) {
282               header_buf[0] = i;
283               header_buf[1] = -1;
284               to = i;
285
286               //start = MPI_Wtime();
287               smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
288
289               /* still need to chop data so that we can use the same non-root code */
290               for (j = 0; j < pipe_length; j++) {
291                 smpi_mpi_send((char *)buf + (j * increment), segment, datatype, to, tag,
292                          comm);
293               }
294
295               //smpi_mpi_send(buf,count,datatype,to,tag,comm);
296               //smpi_mpi_wait(&request,MPI_STATUS_IGNORE);
297
298               //total = MPI_Wtime() - start;
299               //total *= 1000;
300               //printf("SEND TO SINGLE node %d time = %.2f\n",i,total);
301
302
303               already_sent[i] = 1;
304               sent_count++;
305               break;
306             }
307           }
308         }
309
310       }                         /* while loop */
311
312       //total = MPI_Wtime() - start2;
313       //total *= 1000;
314       //printf("Node zero iter = %d time = %.2f\n",iteration,total);
315
316       /* probe before exit in case there are messages to recv */
317       for (i = 1; i < size; i++) {
318         smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
319                    &temp_status_array[i]);
320         if (flag_array[i] == 1)
321           smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
322       }
323     }
324
325     /* rank 0 */
326     /* none root */
327     else {
328
329       /* if root already send a message to this node, don't send one-byte message */
330       smpi_mpi_iprobe(0, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[0], &status);
331
332       /* send 1-byte message to root */
333       if (flag_array[0] == 0)
334         smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
335
336       /* wait for header forward when required */
337       request = smpi_mpi_irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm);
338       smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
339
340       /* search for where it is */
341       int myordering = 0;
342       while (rank != header_buf[myordering]) {
343         myordering++;
344       }
345
346       /* send header when required */
347       if (header_buf[myordering + 1] != -1) {
348         smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
349                  tag, comm);
350       }
351
352       /* receive data */
353
354       if (0 == -1) {
355         //if (header_buf[1] == -1) {
356         request = smpi_mpi_irecv(buf, count, datatype, 0, tag, comm);
357         smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
358         //printf("\t\tnode %d ordering = %d receive data from root\n",rank,myordering);
359       } else {
360         for (i = 0; i < pipe_length; i++) {
361           recv_request_array[i] = smpi_mpi_irecv((char *)buf + (i * increment), segment, datatype, MPI_ANY_SOURCE,
362                     tag, comm);
363         }
364       }
365
366       /* send data */
367       if (header_buf[myordering + 1] != -1) {
368         for (i = 0; i < pipe_length; i++) {
369           smpi_mpi_wait(&recv_request_array[i], MPI_STATUS_IGNORE);
370           send_request_array[i] = smpi_mpi_isend((char *)buf + (i * increment), segment, datatype,
371                     header_buf[myordering + 1], tag, comm);
372         }
373         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
374       }else{
375         for (i = 0; i < pipe_length; i++) 
376           smpi_mpi_wait(&recv_request_array[i], &recv_status_array[i]);
377       }
378
379     }
380
381     free(send_request_array);
382     free(recv_request_array);
383     free(send_status_array);
384     free(recv_status_array);
385   }                             /* end pipeline */
386
387   /* when count is not divisible by block size, use default BCAST for the remainder */
388   if ((remainder != 0) && (count > segment)) {
389     XBT_WARN("MPI_bcast_arrival_nb use default MPI_bcast.");      
390     smpi_mpi_bcast((char *)buf + (pipe_length * increment), remainder, datatype, root, comm);
391   }
392
393   return MPI_SUCCESS;
394 }