Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge pull request #193 from Takishipp/signals
[simgrid.git] / src / smpi / colls / bcast / bcast-arrival-scatter.cpp
1 /* Copyright (c) 2013-2017. The SimGrid Team.
2  * All rights reserved.                                                     */
3
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. */
6
7 #include "../colls_private.h"
8
9 #ifndef BCAST_ARRIVAL_PATTERN_AWARE_HEADER_SIZE
10 #define BCAST_ARRIVAL_PATTERN_AWARE_HEADER_SIZE 128
11 #endif
12
13 #ifndef BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE
14 #define BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE 128
15 #endif
16 namespace simgrid{
17 namespace smpi{
18 /* Non-topology-specific pipelined linear-bcast function */
19 int Coll_bcast_arrival_scatter::bcast(void *buf, int count,
20                                           MPI_Datatype datatype, int root,
21                                           MPI_Comm comm)
22 {
23   int tag = -COLL_TAG_BCAST;//in order to use ANY_TAG, make this one positive
24   int header_tag = 10;
25   MPI_Status status;
26
27   int curr_remainder;
28   int curr_size;
29   int curr_increment;
30   int send_offset;
31   int recv_offset;
32   int send_count;
33   int recv_count;
34
35   MPI_Status temp_status_array[BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE];
36
37   int rank, size;
38   int i, k;
39
40   int sent_count;
41   int header_index;
42   int flag_array[BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE];
43   int already_sent[BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE];
44   int header_buf[BCAST_ARRIVAL_PATTERN_AWARE_HEADER_SIZE];
45   char temp_buf[BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE];
46   int will_send[BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE];
47   int max_node = BCAST_ARRIVAL_PATTERN_AWARE_MAX_NODE;
48   int header_size = BCAST_ARRIVAL_PATTERN_AWARE_HEADER_SIZE;
49
50   MPI_Aint extent;
51   extent = datatype->get_extent();
52
53
54   /* source and destination */
55   int to, from;
56
57   rank = comm->rank();
58   size = comm->size();
59
60   /* message too small */
61   if (count < size) {
62     XBT_WARN("MPI_bcast_arrival_scatter use default MPI_bcast.");
63     Colls::bcast(buf, count, datatype, root, comm);
64     return MPI_SUCCESS;
65   }
66
67
68
69   /* if root is not zero send to rank zero first
70      this can be modified to make it faster by using logical src, dst.
71    */
72   if (root != 0) {
73     if (rank == root) {
74       Request::send(buf, count, datatype, 0, tag - 1, comm);
75     } else if (rank == 0) {
76       Request::recv(buf, count, datatype, root, tag - 1, comm, &status);
77     }
78   }
79
80
81   /* value == 0 means root has not send data (or header) to the node yet */
82   for (i = 0; i < max_node; i++) {
83     already_sent[i] = 0;
84   }
85
86   /* start bcast */
87
88   /* root */
89   if (rank == 0) {
90
91     for (i = 0; i < max_node; i++)
92       will_send[i] = 0;
93
94     sent_count = 0;
95     while (sent_count < (size - 1)) {
96
97       for (k = 0; k < 3; k++) {
98         for (i = 1; i < size; i++) {
99           if ((already_sent[i] == 0) && (will_send[i] == 0)) {
100             Request::iprobe(i, MPI_ANY_TAG, comm, &flag_array[i],
101                        &temp_status_array[i]);
102             if (flag_array[i] == 1) {
103               will_send[i] = 1;
104               Request::recv(&temp_buf[i], 1, MPI_CHAR, i, tag, comm,
105                        &status);
106               i = 0;
107             }
108           }
109         }
110       }
111       header_index = 0;
112
113       /* recv 1-byte message in this round */
114       for (i = 1; i < size; i++) {
115         /* message arrive */
116         if ((will_send[i] == 1) && (already_sent[i] == 0)) {
117           header_buf[header_index] = i;
118           header_index++;
119           sent_count++;
120
121           /* will send in the next step */
122           already_sent[i] = 1;
123         }
124       }
125
126       /*
127          if (header_index != 0) {
128          printf("header index = %d node = ",header_index);
129          for (i=0;i<header_index;i++) {
130          printf("%d ",header_buf[i]);
131          }
132          printf("\n");
133          }
134        */
135
136       /* send header followed by data */
137       if (header_index != 0) {
138         header_buf[header_index] = -1;
139
140         /* send header */
141         for (i = 0; i < header_index; i++) {
142           to = header_buf[i];
143           Request::send(header_buf, header_size, MPI_INT, to, header_tag, comm);
144         }
145
146         curr_remainder = count % header_index;
147         curr_size = (count / header_index);
148         curr_increment = curr_size * extent;
149
150         /* send data */
151
152         for (i = 0; i < header_index; i++) {
153           to = header_buf[i];
154           if ((i == (header_index - 1)) || (curr_size == 0))
155             curr_size += curr_remainder;
156           //printf("Root send to %d index %d\n",to,(i*curr_increment));
157           Request::send((char *) buf + (i * curr_increment), curr_size, datatype, to,
158                    tag, comm);
159         }
160       }
161     }                           /* while (sent_count < size-1) */
162   }
163
164   /* rank 0 */
165   /* none root */
166   else {
167     /* send 1-byte message to root */
168     Request::send(temp_buf, 1, MPI_CHAR, 0, tag, comm);
169
170     /* wait for header forward when required */
171     Request::recv(header_buf, header_size, MPI_INT, 0, header_tag, comm, &status);
172
173     /* search for where it is */
174     int myordering = 0;
175     while (rank != header_buf[myordering]) {
176       myordering++;
177     }
178
179     int total_nodes = 0;
180     while (header_buf[total_nodes] != -1) {
181       total_nodes++;
182     }
183
184     curr_remainder = count % total_nodes;
185     curr_size = (count / total_nodes);
186     curr_increment = curr_size * extent;
187     int recv_size = curr_size;
188
189     /* receive data */
190     if (myordering == (total_nodes - 1))
191       recv_size += curr_remainder;
192     Request::recv((char *) buf + (myordering * curr_increment), recv_size, datatype,
193              0, tag, comm, &status);
194
195     /* at this point all nodes in this set perform all-gather operation */
196
197     to = (myordering == (total_nodes - 1)) ? header_buf[0] : header_buf[myordering + 1];
198     from = (myordering == 0) ? header_buf[total_nodes - 1] : header_buf[myordering - 1];
199
200
201     /* last segment may have a larger size since it also include the remainder */
202     int last_segment_ptr = (total_nodes - 1) * (count / total_nodes) * extent;
203
204
205     /* allgather */
206     for (i = 0; i < total_nodes - 1; i++) {
207       send_offset =
208           ((myordering - i + total_nodes) % total_nodes) * curr_increment;
209       recv_offset =
210           ((myordering - i - 1 + total_nodes) % total_nodes) * curr_increment;
211
212       /* adjust size */
213       if (send_offset != last_segment_ptr)
214         send_count = curr_size;
215       else
216         send_count = curr_size + curr_remainder;
217
218       if (recv_offset != last_segment_ptr)
219         recv_count = curr_size;
220       else
221         recv_count = curr_size + curr_remainder;
222
223       //printf("\t\tnode %d sent_to %d recv_from %d send_size %d recv_size %d\n",rank,to,from,send_count,recv_count);
224       //printf("\tnode %d sent_offset %d send_count %d\n",rank,send_offset,send_count);
225
226
227       Request::sendrecv((char *) buf + send_offset, send_count, datatype, to,
228                    tag + i, (char *) buf + recv_offset, recv_count, datatype,
229                    from, tag + i, comm, &status);
230     }
231   }                             /* non-root */
232
233   return MPI_SUCCESS;
234 }
235
236 }
237 }