1 /* Copyright (c) 2013-2019. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
7 #include "../colls_private.hpp"
8 //#include <star-reduction.c>
10 int reduce_arrival_pattern_aware_segment_size_in_byte = 8192;
13 #define HEADER_SIZE 1024
21 /* Non-topology-specific pipelined linear-reduce function */
22 int Coll_reduce_arrival_pattern_aware::reduce(void *buf, void *rbuf,
24 MPI_Datatype datatype,
28 int rank = comm->rank();
29 int tag = -COLL_TAG_REDUCE;
32 MPI_Request *send_request_array;
33 MPI_Request *recv_request_array;
34 MPI_Status *send_status_array;
35 MPI_Status *recv_status_array;
37 MPI_Status temp_status_array[MAX_NODE];
39 int size = comm->size();
44 int flag_array[MAX_NODE];
45 int already_received[MAX_NODE];
47 int header_buf[HEADER_SIZE];
48 char temp_buf[MAX_NODE];
51 datatype->extent(&lb, &extent);
53 /* source and destination */
56 /* segment is segment size in number of elements (not bytes) */
57 int segment = reduce_arrival_pattern_aware_segment_size_in_byte / extent;
60 int pipe_length = count / segment;
62 /* use for buffer offset for sending and receiving data = segment size in byte */
63 int increment = segment * extent;
65 /* if the input size is not divisible by segment size =>
66 the small remainder will be done with native implementation */
67 int remainder = count % segment;
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_received[i] = 0;
76 tmp_buf = (char *) smpi_get_tmp_sendbuffer(count * extent);
78 Request::sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank,
83 /* when a message is smaller than a block size => no pipeline */
84 if (count <= segment) {
89 while (sent_count < (size - 1)) {
91 for (i = 1; i < size; i++) {
92 if (already_received[i] == 0) {
93 Request::iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
95 simcall_process_sleep(0.0001);
100 /* recv 1-byte message */
101 for (i = 0; i < size; i++) {
105 /* 1-byte message arrive */
106 if ((flag_array[i] == 1) && (already_received[i] == 0)) {
107 Request::recv(temp_buf, 1, MPI_CHAR, i, tag, comm, &status);
108 header_buf[header_index] = i;
113 //printf("root send to %d recv from %d : data = ",to,from);
115 for (i=0;i<=header_index;i++) {
116 printf("%d ",header_buf[i]);
120 /* will receive in the next step */
121 already_received[i] = 1;
125 /* send header followed by receive and reduce data */
126 if (header_index != 0) {
127 header_buf[header_index] = -1;
129 from = header_buf[header_index - 1];
131 Request::send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
132 Request::recv(tmp_buf, count, datatype, from, tag, comm, &status);
133 if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuf, &count, datatype);
142 /* send 1-byte message to root */
143 Request::send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
145 /* wait for header and data, forward when required */
146 Request::recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm,
148 // Request::recv(buf,count,datatype,MPI_ANY_SOURCE,tag,comm,&status);
150 /* search for where it is */
152 while (rank != header_buf[myordering]) {
157 if (header_buf[myordering + 1] != -1) {
158 Request::send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
161 //printf("node %d ordering %d\n",rank,myordering);
163 /* receive, reduce, and forward data */
166 if (myordering == 0) {
167 if (header_buf[myordering + 1] == -1) {
170 to = header_buf[myordering + 1];
172 Request::send(rbuf, count, datatype, to, tag, comm);
175 /* recv, reduce, send */
177 if (header_buf[myordering + 1] == -1) {
180 to = header_buf[myordering + 1];
182 from = header_buf[myordering - 1];
183 Request::recv(tmp_buf, count, datatype, from, tag, comm, &status);
184 if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuf, &count, datatype);
185 Request::send(rbuf, count, datatype, to, tag, comm);
191 // printf("node %d start\n",rank);
194 (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
196 (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
198 (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
200 (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
205 int will_send[MAX_NODE];
206 for (i = 0; i < MAX_NODE; i++)
209 /* loop until all data are received (sent) */
210 while (sent_count < (size - 1)) {
212 for (k = 0; k < 1; k++) {
213 for (i = 1; i < size; i++) {
216 if ((already_received[i] == 0) && (will_send[i] == 0)) {
217 Request::iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
218 &temp_status_array[i]);
219 if (flag_array[i] == 1) {
221 Request::recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
223 //printf("recv from %d\n",i);
228 } /* end of probing */
232 /* recv 1-byte message */
233 for (i = 1; i < size; i++) {
236 /* message arrived in this round (put in the header) */
237 if ((will_send[i] == 1) && (already_received[i] == 0)) {
238 header_buf[header_index] = i;
242 /* will send in the next step */
243 already_received[i] = 1;
247 /* send header followed by data */
248 if (header_index != 0) {
249 header_buf[header_index] = -1;
253 Request::send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm);
255 /* recv data - pipeline */
256 from = header_buf[header_index - 1];
257 for (i = 0; i < pipe_length; i++) {
258 Request::recv(tmp_buf + (i * increment), segment, datatype, from, tag,
260 if(op!=MPI_OP_NULL) op->apply( tmp_buf + (i * increment),
261 (char *)rbuf + (i * increment), &segment, datatype);
264 } /* while loop (sent_count < size-1 ) */
270 /* send 1-byte message to root */
271 Request::send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
274 /* wait for header forward when required */
275 request=Request::irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm);
276 Request::wait(&request, MPI_STATUS_IGNORE);
278 /* search for where it is */
281 while (rank != header_buf[myordering]) {
285 /* send header when required */
286 if (header_buf[myordering + 1] != -1) {
287 Request::send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1],
291 /* (receive, reduce), and send data */
292 if (header_buf[myordering + 1] == -1) {
295 to = header_buf[myordering + 1];
299 if (myordering == 0) {
300 for (i = 0; i < pipe_length; i++) {
301 send_request_array[i]= Request::isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm);
303 Request::waitall((pipe_length), send_request_array, send_status_array);
306 /* receive, reduce, and send */
308 from = header_buf[myordering - 1];
309 for (i = 0; i < pipe_length; i++) {
310 recv_request_array[i]=Request::irecv(tmp_buf + (i * increment), segment, datatype, from, tag, comm);
312 for (i = 0; i < pipe_length; i++) {
313 Request::wait(&recv_request_array[i], MPI_STATUS_IGNORE);
314 if(op!=MPI_OP_NULL) op->apply( tmp_buf + (i * increment), (char *)rbuf + (i * increment),
316 send_request_array[i]=Request::isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm);
318 Request::waitall((pipe_length), send_request_array, send_status_array);
325 free(send_request_array);
326 free(recv_request_array);
327 free(send_status_array);
328 free(recv_status_array);
330 //printf("node %d done\n",rank);
334 /* if root is not zero send root after finished
335 this can be modified to make it faster by using logical src, dst.
339 Request::send(rbuf, count, datatype, root, tag, comm);
340 } else if (rank == root) {
341 Request::recv(rbuf, count, datatype, 0, tag, comm, &status);
346 /* when count is not divisible by block size, use default BCAST for the remainder */
347 if ((remainder != 0) && (count > segment)) {
348 Coll_reduce_default::reduce((char*)buf + (pipe_length * increment), (char*)rbuf + (pipe_length * increment),
349 remainder, datatype, op, root, comm);
352 smpi_free_tmp_buffer(tmp_buf);