Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
kill all trailling whitespaces
[simgrid.git] / src / smpi / colls / barrier / barrier-ompi.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 /*
8  * Copyright (c) 2004-2005 The Trustees of Indiana University and Indiana
9  *                         University Research and Technology
10  *                         Corporation.  All rights reserved.
11  * Copyright (c) 2004-2006 The University of Tennessee and The University
12  *                         of Tennessee Research Foundation.  All rights
13  *                         reserved.
14  * Copyright (c) 2004-2005 High Performance Computing Center Stuttgart,
15  *                         University of Stuttgart.  All rights reserved.
16  * Copyright (c) 2004-2005 The Regents of the University of California.
17  *                         All rights reserved.
18  * Copyright (c) 2008      Sun Microsystems, Inc.  All rights reserved.
19  *
20  * Additional copyrights may follow
21  */
22
23 #include "../colls_private.h"
24 #include "../coll_tuned_topo.h"
25
26
27 /*
28  * Barrier is ment to be a synchronous operation, as some BTLs can mark
29  * a request done before its passed to the NIC and progress might not be made
30  * elsewhere we cannot allow a process to exit the barrier until its last
31  * [round of] sends are completed.
32  *
33  * It is last round of sends rather than 'last' individual send as each pair of
34  * peers can use different channels/devices/btls and the receiver of one of
35  * these sends might be forced to wait as the sender
36  * leaves the collective and does not make progress until the next mpi call
37  *
38  */
39
40 /*
41  * Simple double ring version of barrier
42  *
43  * synchronous gurantee made by last ring of sends are synchronous
44  *
45  */
46 namespace simgrid{
47 namespace smpi{
48 int Coll_barrier_ompi_doublering::barrier(MPI_Comm comm)
49 {
50     int rank, size;
51     int left, right;
52
53
54     rank = comm->rank();
55     size = comm->size();
56
57     XBT_DEBUG("ompi_coll_tuned_barrier_ompi_doublering rank %d", rank);
58
59     left = ((rank-1+size)%size);
60     right = ((rank+1)%size);
61
62     if (rank > 0) { /* receive message from the left */
63         Request::recv((void*)NULL, 0, MPI_BYTE, left,
64                                 COLL_TAG_BARRIER, comm,
65                                 MPI_STATUS_IGNORE);
66     }
67
68     /* Send message to the right */
69     Request::send((void*)NULL, 0, MPI_BYTE, right,
70                             COLL_TAG_BARRIER,
71                              comm);
72
73     /* root needs to receive from the last node */
74     if (rank == 0) {
75         Request::recv((void*)NULL, 0, MPI_BYTE, left,
76                                 COLL_TAG_BARRIER, comm,
77                                 MPI_STATUS_IGNORE);
78     }
79
80     /* Allow nodes to exit */
81     if (rank > 0) { /* post Receive from left */
82         Request::recv((void*)NULL, 0, MPI_BYTE, left,
83                                 COLL_TAG_BARRIER, comm,
84                                 MPI_STATUS_IGNORE);
85     }
86
87     /* send message to the right one */
88     Request::send((void*)NULL, 0, MPI_BYTE, right,
89                             COLL_TAG_BARRIER,
90                              comm);
91
92     /* rank 0 post receive from the last node */
93     if (rank == 0) {
94         Request::recv((void*)NULL, 0, MPI_BYTE, left,
95                                 COLL_TAG_BARRIER, comm,
96                                 MPI_STATUS_IGNORE);
97     }
98
99     return MPI_SUCCESS;
100
101 }
102
103
104 /*
105  * To make synchronous, uses sync sends and sync sendrecvs
106  */
107
108 int Coll_barrier_ompi_recursivedoubling::barrier(MPI_Comm comm)
109 {
110     int rank, size, adjsize;
111     int mask, remote;
112
113     rank = comm->rank();
114     size = comm->size();
115     XBT_DEBUG(
116                  "ompi_coll_tuned_barrier_ompi_recursivedoubling rank %d",
117                  rank);
118
119     /* do nearest power of 2 less than size calc */
120     for( adjsize = 1; adjsize <= size; adjsize <<= 1 );
121     adjsize >>= 1;
122
123     /* if size is not exact power of two, perform an extra step */
124     if (adjsize != size) {
125         if (rank >= adjsize) {
126             /* send message to lower ranked node */
127             remote = rank - adjsize;
128             Request::sendrecv(NULL, 0, MPI_BYTE, remote,
129                                                   COLL_TAG_BARRIER,
130                                                   NULL, 0, MPI_BYTE, remote,
131                                                   COLL_TAG_BARRIER,
132                                                   comm, MPI_STATUS_IGNORE);
133
134         } else if (rank < (size - adjsize)) {
135
136             /* receive message from high level rank */
137             Request::recv((void*)NULL, 0, MPI_BYTE, rank+adjsize,
138                                     COLL_TAG_BARRIER, comm,
139                                     MPI_STATUS_IGNORE);
140
141         }
142     }
143
144     /* exchange messages */
145     if ( rank < adjsize ) {
146         mask = 0x1;
147         while ( mask < adjsize ) {
148             remote = rank ^ mask;
149             mask <<= 1;
150             if (remote >= adjsize) continue;
151
152             /* post receive from the remote node */
153             Request::sendrecv(NULL, 0, MPI_BYTE, remote,
154                                                   COLL_TAG_BARRIER,
155                                                   NULL, 0, MPI_BYTE, remote,
156                                                   COLL_TAG_BARRIER,
157                                                   comm, MPI_STATUS_IGNORE);
158         }
159     }
160
161     /* non-power of 2 case */
162     if (adjsize != size) {
163         if (rank < (size - adjsize)) {
164             /* send enter message to higher ranked node */
165             remote = rank + adjsize;
166             Request::send((void*)NULL, 0, MPI_BYTE, remote,
167                                     COLL_TAG_BARRIER,
168                                      comm);
169
170         }
171     }
172
173     return MPI_SUCCESS;
174
175 }
176
177
178 /*
179  * To make synchronous, uses sync sends and sync sendrecvs
180  */
181
182 int Coll_barrier_ompi_bruck::barrier(MPI_Comm comm)
183 {
184     int rank, size;
185     int distance, to, from;
186
187     rank = comm->rank();
188     size = comm->size();
189     XBT_DEBUG(
190                  "ompi_coll_tuned_barrier_ompi_bruck rank %d", rank);
191
192     /* exchange data with rank-2^k and rank+2^k */
193     for (distance = 1; distance < size; distance <<= 1) {
194         from = (rank + size - distance) % size;
195         to   = (rank + distance) % size;
196
197         /* send message to lower ranked node */
198         Request::sendrecv(NULL, 0, MPI_BYTE, to,
199                                               COLL_TAG_BARRIER,
200                                               NULL, 0, MPI_BYTE, from,
201                                               COLL_TAG_BARRIER,
202                                               comm, MPI_STATUS_IGNORE);
203     }
204
205     return MPI_SUCCESS;
206
207 }
208
209
210 /*
211  * To make synchronous, uses sync sends and sync sendrecvs
212  */
213 /* special case for two processes */
214 int Coll_barrier_ompi_two_procs::barrier(MPI_Comm comm)
215 {
216     int remote;
217
218     remote = comm->rank();
219     XBT_DEBUG(
220                  "ompi_coll_tuned_barrier_ompi_two_procs rank %d", remote);
221     remote = (remote + 1) & 0x1;
222
223     Request::sendrecv(NULL, 0, MPI_BYTE, remote,
224                                           COLL_TAG_BARRIER,
225                                           NULL, 0, MPI_BYTE, remote,
226                                           COLL_TAG_BARRIER,
227                                           comm, MPI_STATUS_IGNORE);
228     return (MPI_SUCCESS);
229 }
230
231
232 /*
233  * Linear functions are copied from the BASIC coll module
234  * they do not segment the message and are simple implementations
235  * but for some small number of nodes and/or small data sizes they
236  * are just as fast as tuned/tree based segmenting operations
237  * and as such may be selected by the decision functions
238  * These are copied into this module due to the way we select modules
239  * in V1. i.e. in V2 we will handle this differently and so will not
240  * have to duplicate code.
241  * GEF Oct05 after asking Jeff.
242  */
243
244 /* copied function (with appropriate renaming) starts here */
245
246 int Coll_barrier_ompi_basic_linear::barrier(MPI_Comm comm)
247 {
248     int i;
249     int size = comm->size();
250     int rank = comm->rank();
251
252     /* All non-root send & receive zero-length message. */
253
254     if (rank > 0) {
255         Request::send (NULL, 0, MPI_BYTE, 0,
256                                  COLL_TAG_BARRIER,
257                                   comm);
258
259         Request::recv (NULL, 0, MPI_BYTE, 0,
260                                  COLL_TAG_BARRIER,
261                                  comm, MPI_STATUS_IGNORE);
262     }
263
264     /* The root collects and broadcasts the messages. */
265
266     else {
267         MPI_Request* requests;
268
269         requests = (MPI_Request*)malloc( size * sizeof(MPI_Request) );
270         for (i = 1; i < size; ++i) {
271             requests[i] = Request::irecv(NULL, 0, MPI_BYTE, MPI_ANY_SOURCE,
272                                      COLL_TAG_BARRIER, comm
273                                      );
274         }
275         Request::waitall( size-1, requests+1, MPI_STATUSES_IGNORE );
276
277         for (i = 1; i < size; ++i) {
278             requests[i] = Request::isend(NULL, 0, MPI_BYTE, i,
279                                      COLL_TAG_BARRIER,
280                                       comm
281                                      );
282         }
283         Request::waitall( size-1, requests+1, MPI_STATUSES_IGNORE );
284         free( requests );
285     }
286
287     /* All done */
288
289     return MPI_SUCCESS;
290
291 }
292 /* copied function (with appropriate renaming) ends here */
293
294 /*
295  * Another recursive doubling type algorithm, but in this case
296  * we go up the tree and back down the tree.
297  */
298 int Coll_barrier_ompi_tree::barrier(MPI_Comm comm)
299 {
300     int rank, size, depth;
301     int jump, partner;
302
303     rank = comm->rank();
304     size = comm->size();
305     XBT_DEBUG(
306                  "ompi_coll_tuned_barrier_ompi_tree %d",
307                  rank);
308
309     /* Find the nearest power of 2 of the communicator size. */
310     for(depth = 1; depth < size; depth <<= 1 );
311
312     for (jump=1; jump<depth; jump<<=1) {
313         partner = rank ^ jump;
314         if (!(partner & (jump-1)) && partner < size) {
315             if (partner > rank) {
316                 Request::recv (NULL, 0, MPI_BYTE, partner,
317                                          COLL_TAG_BARRIER, comm,
318                                          MPI_STATUS_IGNORE);
319             } else if (partner < rank) {
320                 Request::send (NULL, 0, MPI_BYTE, partner,
321                                          COLL_TAG_BARRIER,
322                                           comm);
323             }
324         }
325     }
326
327     depth>>=1;
328     for (jump = depth; jump>0; jump>>=1) {
329         partner = rank ^ jump;
330         if (!(partner & (jump-1)) && partner < size) {
331             if (partner > rank) {
332                 Request::send (NULL, 0, MPI_BYTE, partner,
333                                          COLL_TAG_BARRIER,
334                                           comm);
335             } else if (partner < rank) {
336                 Request::recv (NULL, 0, MPI_BYTE, partner,
337                                          COLL_TAG_BARRIER, comm,
338                                          MPI_STATUS_IGNORE);
339             }
340         }
341     }
342
343     return MPI_SUCCESS;
344 }
345
346 }
347 }