1 #include "colls_private.h"
3 static int bcast_NTSL_segment_size_in_byte = 8192;
5 #define HEADER_SIZE 1024
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,
13 int tag = -COLL_TAG_BCAST;
16 MPI_Request *send_request_array;
17 MPI_Request *recv_request_array;
18 MPI_Status *send_status_array;
19 MPI_Status *recv_status_array;
21 MPI_Status temp_status_array[MAX_NODE];
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];
35 extent = smpi_datatype_get_extent(datatype);
42 rank = smpi_comm_rank(comm);
43 size = smpi_comm_size(comm);
46 /* segment is segment size in number of elements (not bytes) */
47 int segment = bcast_NTSL_segment_size_in_byte / extent;
50 int pipe_length = count / segment;
52 /* use for buffer offset for sending and receiving data = segment size in byte */
53 int increment = segment * extent;
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;
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.
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);
70 /* value == 0 means root has not send data (or header) to the node yet */
71 for (i = 0; i < MAX_NODE; i++) {
76 /* when a message is smaller than a block size => no pipeline */
77 if (count <= segment) {
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],
88 /* recv 1-byte message */
89 for (i = 1; i < size; i++) {
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;
98 /* will send in the next step */
103 /* send header followed by data */
104 if (header_index != 0) {
105 header_buf[header_index] = -1;
107 smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
108 smpi_mpi_send(buf, count, datatype, to, tag, comm);
111 /* randomly MPI_Send to one */
113 /* search for the first node that never received data before */
114 for (i = 1; i < size; i++) {
115 if (already_sent[i] == 0) {
118 smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, i, tag, comm);
119 smpi_mpi_send(buf, count, datatype, i, tag, comm);
134 /* send 1-byte message to root */
135 smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
137 /* wait for header and data, forward when required */
138 smpi_mpi_recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
140 smpi_mpi_recv(buf, count, datatype, MPI_ANY_SOURCE, tag, comm, &status);
142 /* search for where it is */
144 while (rank != header_buf[myordering]) {
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],
152 smpi_mpi_send(buf, count, datatype, header_buf[myordering + 1], tag, comm);
159 (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
161 (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
163 (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
165 (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
168 //double start2 = MPI_Wtime();
171 while (sent_count < (size - 1)) {
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]);
178 //total = MPI_Wtime() - start;
180 //printf("Iprobe time = %.2f\n",total);
184 /* recv 1-byte message */
185 for (i = 1; i < size; i++) {
187 if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
188 smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
190 header_buf[header_index] = i;
194 /* will send in the next step */
198 //total = MPI_Wtime() - start;
200 //printf("Recv 1-byte time = %.2f\n",total);
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]);
212 /* send header followed by data */
213 if (header_index != 0) {
214 header_buf[header_index] = -1;
217 //start = MPI_Wtime();
220 smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
222 //total = MPI_Wtime() - start;
224 //printf("\tSend header to %d time = %.2f\n",to,total);
226 //start = MPI_Wtime();
228 /* send data - non-pipeline case */
231 //if (header_index == 1) {
232 smpi_mpi_send(buf, count, datatype, to, tag, comm);
236 /* send data - pipeline */
238 for (i = 0; i < pipe_length; i++) {
239 smpi_mpi_send((char *)buf + (i * increment), segment, datatype, to, tag, comm);
241 //smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
243 //total = MPI_Wtime() - start;
245 //printf("\tSend data to %d time = %.2f\n",to,total);
251 /* randomly MPI_Send to one node */
253 /* search for the first node that never received data before */
254 for (i = 1; i < size; i++) {
255 if (already_sent[i] == 0) {
260 //start = MPI_Wtime();
261 smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
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,
269 //smpi_mpi_send(buf,count,datatype,to,tag,comm);
270 //smpi_mpi_wait(&request,MPI_STATUS_IGNORE);
272 //total = MPI_Wtime() - start;
274 //printf("SEND TO SINGLE node %d time = %.2f\n",i,total);
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,
290 //total = MPI_Wtime() - start2;
292 //printf("Node zero iter = %d time = %.2f\n",iteration,total);
298 /* send 1-byte message to root */
299 smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
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);
305 /* search for where it is */
307 while (rank != header_buf[myordering]) {
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],
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);
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,
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);
338 smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
340 smpi_mpi_waitall(pipe_length, recv_request_array, recv_status_array);
345 free(send_request_array);
346 free(recv_request_array);
347 free(send_status_array);
348 free(recv_status_array);
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);