Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add tesh files to test all new collectives
[simgrid.git] / src / smpi / colls / bcast-arrival-nb.c
1 #include "colls.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 = 50;
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
31   int header_buf[HEADER_SIZE];
32   char temp_buf[MAX_NODE];
33
34   MPI_Aint extent;
35   MPI_Type_extent(datatype, &extent);
36
37   /* destination */
38   int to;
39
40
41
42   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
43   MPI_Comm_size(MPI_COMM_WORLD, &size);
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       MPI_Send(buf, count, datatype, 0, tag, comm);
65     } else if (rank == 0) {
66       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   }
74   //  printf("YYY\n");
75
76   /* when a message is smaller than a block size => no pipeline */
77   if (count <= segment) {
78     if (rank == 0) {
79       sent_count = 0;
80
81       while (sent_count < (size - 1)) {
82
83         //      for (j=0;j<1000;j++) {
84         for (i = 1; i < size; i++) {
85           if (already_sent[i] == 0)
86             MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
87                        MPI_STATUSES_IGNORE);
88         }
89         //}
90
91         header_index = 0;
92         /* recv 1-byte message */
93         for (i = 1; i < size; i++) {
94
95           /* message arrive */
96           if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
97             MPI_Recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
98             header_buf[header_index] = i;
99             header_index++;
100             sent_count++;
101
102             /* will send in the next step */
103             already_sent[i] = 1;
104           }
105         }
106
107         /* send header followed by data */
108         if (header_index != 0) {
109           header_buf[header_index] = -1;
110           to = header_buf[0];
111           MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
112           MPI_Send(buf, count, datatype, to, tag, comm);
113         }
114
115         /* randomly MPI_Send to one */
116         else {
117           /* search for the first node that never received data before */
118           for (i = 1; i < size; i++) {
119             if (already_sent[i] == 0) {
120               header_buf[0] = i;
121               header_buf[1] = -1;
122               MPI_Send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
123               MPI_Send(buf, count, datatype, i, tag, comm);
124               already_sent[i] = 1;
125               sent_count++;
126               break;
127             }
128           }
129         }
130
131
132       }                         /* while loop */
133     }
134
135     /* non-root */
136     else {
137
138       /* send 1-byte message to root */
139       MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
140
141       /* wait for header and data, forward when required */
142       MPI_Recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
143                &status);
144       MPI_Recv(buf, count, datatype, MPI_ANY_SOURCE, tag, comm, &status);
145
146       /* search for where it is */
147       int myordering = 0;
148       while (rank != header_buf[myordering]) {
149         myordering++;
150       }
151
152       /* send header followed by data */
153       if (header_buf[myordering + 1] != -1) {
154         MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
155                  tag, comm);
156         MPI_Send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
157       }
158     }
159   }
160   /* pipeline bcast */
161   else {
162     send_request_array =
163         (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
164     recv_request_array =
165         (MPI_Request *) malloc((size + pipe_length) * sizeof(MPI_Request));
166     send_status_array =
167         (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
168     recv_status_array =
169         (MPI_Status *) malloc((size + pipe_length) * sizeof(MPI_Status));
170
171     if (rank == 0) {
172       sent_count = 0;
173       int iteration = 0;
174
175       int will_send[1000];
176       for (i = 0; i < 1000; i++)
177         will_send[i] = 0;
178       while (sent_count < (size - 1)) {
179         iteration++;
180         //start = MPI_Wtime();
181
182         int k;
183         for (k = 0; k < 3; k++) {
184           for (i = 1; i < size; i++) {
185             if ((already_sent[i] == 0) && (will_send[i] == 0)) {
186               MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
187                          &temp_status_array[i]);
188               if (flag_array[i] == 1) {
189                 will_send[i] = 1;
190                 MPI_Recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD,
191                          &status);
192                 i = 1;
193               }
194             }
195           }
196         }
197
198         //total = MPI_Wtime() - start;
199         //total *= 1000;
200         //printf("Iprobe time = %.2f\n",total);
201         header_index = 0;
202
203         //start = MPI_Wtime();
204         /* recv 1-byte message */
205         for (i = 1; i < size; i++) {
206           /* message arrive */
207           if ((will_send[i] == 1) && (already_sent[i] == 0)) {
208             header_buf[header_index] = i;
209             header_index++;
210             sent_count++;
211
212             /* will send in the next step */
213             already_sent[i] = 1;
214           }
215         }
216         //printf("sent_count = %d\n",sent_count);
217
218
219         //total = MPI_Wtime() - start;
220         //total *= 1000;
221         //printf("Recv 1-byte time = %.2f\n",total);
222
223         /*      
224            if (header_index != 0) {
225            printf("header index = %d node = ",header_index);
226            for (i=0;i<header_index;i++) {
227            printf("%d ",header_buf[i]);
228            }
229            printf("\n");
230            }
231          */
232
233         /* send header followed by data */
234         if (header_index != 0) {
235           header_buf[header_index] = -1;
236           to = header_buf[0];
237
238           //start = MPI_Wtime();
239
240           /* send header */
241           MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
242
243           //total = MPI_Wtime() - start;
244           //total *= 1000;
245           //printf("\tSend header to %d time = %.2f\n",to,total);
246
247           //start = MPI_Wtime();
248
249           /* send data - non-pipeline case */
250
251           if (0 == 1) {
252             //if (header_index == 1) {
253             MPI_Send(buf, count, datatype, to, tag, comm);
254           }
255
256
257           /* send data - pipeline */
258           else {
259             for (i = 0; i < pipe_length; i++) {
260               MPI_Send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
261             }
262             //MPI_Waitall((pipe_length), send_request_array, send_status_array);
263           }
264           //total = MPI_Wtime() - start;
265           //total *= 1000;
266           //printf("\tSend data to %d time = %.2f\n",to,total);
267
268         }
269
270
271
272         /* randomly MPI_Send to one node */
273         else {
274           /* search for the first node that never received data before */
275           for (i = 1; i < size; i++) {
276             if (already_sent[i] == 0) {
277               header_buf[0] = i;
278               header_buf[1] = -1;
279               to = i;
280
281               //start = MPI_Wtime();
282               MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
283
284               /* still need to chop data so that we can use the same non-root code */
285               for (j = 0; j < pipe_length; j++) {
286                 MPI_Send((char *)buf + (j * increment), segment, datatype, to, tag,
287                          comm);
288               }
289
290               //MPI_Send(buf,count,datatype,to,tag,comm);
291               //MPI_Wait(&request,MPI_STATUS_IGNORE);
292
293               //total = MPI_Wtime() - start;
294               //total *= 1000;
295               //printf("SEND TO SINGLE node %d time = %.2f\n",i,total);
296
297
298               already_sent[i] = 1;
299               sent_count++;
300               break;
301             }
302           }
303         }
304
305       }                         /* while loop */
306
307       //total = MPI_Wtime() - start2;
308       //total *= 1000;
309       //printf("Node zero iter = %d time = %.2f\n",iteration,total);
310
311       /* probe before exit in case there are messages to recv */
312       for (i = 1; i < size; i++) {
313         MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i],
314                    &temp_status_array[i]);
315         if (flag_array[i] == 1)
316           MPI_Recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
317       }
318     }
319
320     /* rank 0 */
321     /* none root */
322     else {
323
324       /* if root already send a message to this node, don't send one-byte message */
325       MPI_Iprobe(0, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[0], &status);
326
327       /* send 1-byte message to root */
328       if (flag_array[0] == 0)
329         MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
330
331       /* wait for header forward when required */
332       MPI_Irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
333                 &request);
334       MPI_Wait(&request, MPI_STATUS_IGNORE);
335
336       /* search for where it is */
337       int myordering = 0;
338       while (rank != header_buf[myordering]) {
339         myordering++;
340       }
341
342       /* send header when required */
343       if (header_buf[myordering + 1] != -1) {
344         MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
345                  tag, comm);
346       }
347
348       /* receive data */
349
350       if (0 == -1) {
351         //if (header_buf[1] == -1) {
352         MPI_Irecv(buf, count, datatype, 0, tag, comm, &request);
353         MPI_Wait(&request, MPI_STATUS_IGNORE);
354         //printf("\t\tnode %d ordering = %d receive data from root\n",rank,myordering);
355       } else {
356         for (i = 0; i < pipe_length; i++) {
357           MPI_Irecv((char *)buf + (i * increment), segment, datatype, MPI_ANY_SOURCE,
358                     tag, comm, &recv_request_array[i]);
359         }
360       }
361
362       /* send data */
363       if (header_buf[myordering + 1] != -1) {
364         for (i = 0; i < pipe_length; i++) {
365           MPI_Wait(&recv_request_array[i], MPI_STATUS_IGNORE);
366           MPI_Isend((char *)buf + (i * increment), segment, datatype,
367                     header_buf[myordering + 1], tag, comm,
368                     &send_request_array[i]);
369         }
370         MPI_Waitall((pipe_length), send_request_array, send_status_array);
371       }
372
373     }
374
375     free(send_request_array);
376     free(recv_request_array);
377     free(send_status_array);
378     free(recv_status_array);
379   }                             /* end pipeline */
380
381   /* when count is not divisible by block size, use default BCAST for the remainder */
382   if ((remainder != 0) && (count > segment)) {
383     MPI_Bcast((char *)buf + (pipe_length * increment), remainder, datatype, root, comm);
384   }
385
386   return MPI_SUCCESS;
387 }