Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
use tuned algo here
[simgrid.git] / src / smpi / colls / bcast-arrival-pattern-aware.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_pattern_aware(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(comm);
43   size = smpi_comm_size(comm);
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
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         for (i = 1; i < size; i++) {
83           smpi_mpi_iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
84                      MPI_STATUSES_IGNORE);
85         }
86
87         header_index = 0;
88         /* recv 1-byte message */
89         for (i = 1; i < size; i++) {
90
91           /* message arrive */
92           if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
93             smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, comm, &status);
94             header_buf[header_index] = i;
95             header_index++;
96             sent_count++;
97
98             /* will send in the next step */
99             already_sent[i] = 1;
100           }
101         }
102
103         /* send header followed by data */
104         if (header_index != 0) {
105           header_buf[header_index] = -1;
106           to = header_buf[0];
107           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
108           smpi_mpi_send(buf, count, datatype, to, tag, comm);
109         }
110
111         /* randomly MPI_Send to one */
112         else {
113           /* search for the first node that never received data before */
114           for (i = 1; i < size; i++) {
115             if (already_sent[i] == 0) {
116               header_buf[0] = i;
117               header_buf[1] = -1;
118               smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
119               smpi_mpi_send(buf, count, datatype, i, tag, comm);
120               already_sent[i] = 1;
121               sent_count++;
122               break;
123             }
124           }
125         }
126
127
128       }                         /* while loop */
129     }
130
131     /* non-root */
132     else {
133
134       /* send 1-byte message to root */
135       smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
136
137       /* wait for header and data, forward when required */
138       smpi_mpi_recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
139                &status);
140       smpi_mpi_recv(buf, count, datatype, MPI_ANY_SOURCE, tag, comm, &status);
141
142       /* search for where it is */
143       int myordering = 0;
144       while (rank != header_buf[myordering]) {
145         myordering++;
146       }
147
148       /* send header followed by data */
149       if (header_buf[myordering + 1] != -1) {
150         smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
151                  tag, comm);
152         smpi_mpi_send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
153       }
154     }
155   }
156   /* pipeline bcast */
157   else {
158     send_request_array =
159         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
160     recv_request_array =
161         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
162     send_status_array =
163         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
164     recv_status_array =
165         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
166
167     if (rank == 0) {
168       //double start2 = MPI_Wtime();
169       sent_count = 0;
170       //int iteration = 0;
171       while (sent_count < (size - 1)) {
172         //iteration++;
173         //start = MPI_Wtime();
174         for (i = 1; i < size; i++) {
175           smpi_mpi_iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
176                      &temp_status_array[i]);
177         }
178         //total = MPI_Wtime() - start;
179         //total *= 1000;
180         //printf("Iprobe time = %.2f\n",total);
181         header_index = 0;
182
183         MPI_Wtime();
184         /* recv 1-byte message */
185         for (i = 1; i < size; i++) {
186           /* message arrive */
187           if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
188             smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
189                      &status);
190             header_buf[header_index] = i;
191             header_index++;
192             sent_count++;
193
194             /* will send in the next step */
195             already_sent[i] = 1;
196           }
197         }
198         //total = MPI_Wtime() - start;
199         //total *= 1000;
200         //printf("Recv 1-byte time = %.2f\n",total);
201
202         /*
203            if (header_index != 0) {
204            printf("header index = %d node = ",header_index);
205            for (i=0;i<header_index;i++) {
206            printf("%d ",header_buf[i]);
207            }
208            printf("\n");
209            }
210          */
211
212         /* send header followed by data */
213         if (header_index != 0) {
214           header_buf[header_index] = -1;
215           to = header_buf[0];
216
217           //start = MPI_Wtime();
218
219           /* send header */
220           smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
221
222           //total = MPI_Wtime() - start;
223           //total *= 1000;
224           //printf("\tSend header to %d time = %.2f\n",to,total);
225
226           //start = MPI_Wtime();
227
228           /* send data - non-pipeline case */
229
230           if (0 == 1) {
231             //if (header_index == 1) {
232             smpi_mpi_send(buf, count, datatype, to, tag, comm);
233           }
234
235
236           /* send data - pipeline */
237           else {
238             for (i = 0; i < pipe_length; i++) {
239               smpi_mpi_send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
240             }
241             //smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
242           }
243           //total = MPI_Wtime() - start;
244           //total *= 1000;
245           //printf("\tSend data to %d time = %.2f\n",to,total);
246
247         }
248
249
250
251         /* randomly MPI_Send to one node */
252         else {
253           /* search for the first node that never received data before */
254           for (i = 1; i < size; i++) {
255             if (already_sent[i] == 0) {
256               header_buf[0] = i;
257               header_buf[1] = -1;
258               to = i;
259
260               //start = MPI_Wtime();
261               smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
262
263               /* still need to chop data so that we can use the same non-root code */
264               for (j = 0; j < pipe_length; j++) {
265                 smpi_mpi_send((char *)buf + (j * increment), segment, datatype, to, tag,
266                          comm);
267               }
268
269               //smpi_mpi_send(buf,count,datatype,to,tag,comm);
270               //smpi_mpi_wait(&request,MPI_STATUS_IGNORE);
271
272               //total = MPI_Wtime() - start;
273               //total *= 1000;
274               //printf("SEND TO SINGLE node %d time = %.2f\n",i,total);
275
276
277               already_sent[i] = 1;
278               to_clean[i]=1;
279               sent_count++;
280               break;
281             }
282           }
283         }
284
285       }                         /* while loop */
286
287       for(i=0; i<size; i++)
288         if(to_clean[i]!=0)smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
289                      &status);
290       //total = MPI_Wtime() - start2;
291       //total *= 1000;
292       //printf("Node zero iter = %d time = %.2f\n",iteration,total);
293     }
294
295     /* rank 0 */
296     /* none root */
297     else {
298       /* send 1-byte message to root */
299       smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
300
301       /* wait for header forward when required */
302       request = smpi_mpi_irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm);
303       smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
304
305       /* search for where it is */
306       int myordering = 0;
307       while (rank != header_buf[myordering]) {
308         myordering++;
309       }
310
311       /* send header when required */
312       if (header_buf[myordering + 1] != -1) {
313         smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
314                  tag, comm);
315       }
316
317       /* receive data */
318
319       if (0 == -1) {
320         //if (header_buf[1] == -1) {
321         request = smpi_mpi_irecv(buf, count, datatype, 0, tag, comm);
322         smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
323         //printf("\t\tnode %d ordering = %d receive data from root\n",rank,myordering);
324       } else {
325         for (i = 0; i < pipe_length; i++) {
326           recv_request_array[i] = smpi_mpi_irecv((char *)buf + (i * increment), segment, datatype, MPI_ANY_SOURCE,
327                                                  tag, comm);
328         }
329       }
330
331       /* send data */
332       if (header_buf[myordering + 1] != -1) {
333         for (i = 0; i < pipe_length; i++) {
334           smpi_mpi_wait(&recv_request_array[i], MPI_STATUS_IGNORE);
335           send_request_array[i] = smpi_mpi_isend((char *)buf + (i * increment), segment, datatype,
336                     header_buf[myordering + 1], tag, comm);
337         }
338         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
339       }else{
340           smpi_mpi_waitall(pipe_length, recv_request_array, recv_status_array);
341           }
342     
343     }
344
345     free(send_request_array);
346     free(recv_request_array);
347     free(send_status_array);
348     free(recv_status_array);
349   }                             /* end pipeline */
350
351   /* when count is not divisible by block size, use default BCAST for the remainder */
352   if ((remainder != 0) && (count > segment)) {
353     XBT_WARN("MPI_bcast_arrival_pattern_aware use default MPI_bcast.");   
354     smpi_mpi_bcast((char *)buf + (pipe_length * increment), remainder, datatype, root, comm);
355   }
356
357   return MPI_SUCCESS;
358 }