Logo AND Algorithmique Numérique Distribuée

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