Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Sort include lists according to clang-format.
[simgrid.git] / src / smpi / colls / gather / gather-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-2009 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  *
19  * Additional copyrights may follow
20  */
21
22 #include "../coll_tuned_topo.hpp"
23 #include "../colls_private.hpp"
24
25 namespace simgrid{
26 namespace smpi{
27
28 int Coll_gather_ompi_binomial::gather(void* sbuf, int scount, MPI_Datatype sdtype, void* rbuf, int rcount,
29                                       MPI_Datatype rdtype, int root, MPI_Comm comm)
30 {
31     int line = -1;
32     int i;
33     int rank;
34     int vrank;
35     int size;
36     int total_recv = 0;
37     char *ptmp     = NULL;
38     char *tempbuf  = NULL;
39     int err;
40     ompi_coll_tree_t* bmtree;
41     MPI_Status status;
42     MPI_Aint sextent, slb, strue_lb, strue_extent;
43     MPI_Aint rextent, rlb, rtrue_lb, rtrue_extent;
44
45
46     size = comm->size();
47     rank = comm->rank();
48
49     XBT_DEBUG("smpi_coll_tuned_gather_ompi_binomial rank %d", rank);
50
51     /* create the binomial tree */
52    // COLL_TUNED_UPDATE_IN_ORDER_BMTREE( comm, tuned_module, root );
53     bmtree = ompi_coll_tuned_topo_build_in_order_bmtree(comm, root);
54     // data->cached_in_order_bmtree;
55
56     sdtype->extent(&slb, &sextent);
57     sdtype->extent(&strue_lb, &strue_extent);
58
59     vrank = (rank - root + size) % size;
60
61     if (rank == root) {
62         rdtype->extent(&rlb, &rextent);
63         rdtype->extent(&rtrue_lb, &rtrue_extent);
64         if (0 == root) {
65           /* root on 0, just use the recv buffer */
66           ptmp = (char*)rbuf;
67           if (sbuf != MPI_IN_PLACE) {
68             err = Datatype::copy(sbuf, scount, sdtype, ptmp, rcount, rdtype);
69             if (MPI_SUCCESS != err) {
70               line = __LINE__;
71               goto err_hndl;
72             }
73           }
74         } else {
75           /* root is not on 0, allocate temp buffer for recv,
76            * rotate data at the end */
77           tempbuf = (char*)smpi_get_tmp_recvbuffer(rtrue_extent + (rcount * size - 1) * rextent);
78           if (NULL == tempbuf) {
79             err  = MPI_ERR_OTHER;
80             line = __LINE__;
81             goto err_hndl;
82           }
83
84           ptmp = tempbuf - rlb;
85           if (sbuf != MPI_IN_PLACE) {
86             /* copy from sbuf to temp buffer */
87             err = Datatype::copy(sbuf, scount, sdtype, ptmp, rcount, rdtype);
88             if (MPI_SUCCESS != err) {
89               line = __LINE__;
90               goto err_hndl;
91             }
92           } else {
93             /* copy from rbuf to temp buffer  */
94             err = Datatype::copy((char*)rbuf + rank * rextent * rcount, rcount, rdtype, ptmp, rcount, rdtype);
95             if (MPI_SUCCESS != err) {
96               line = __LINE__;
97               goto err_hndl;
98             }
99           }
100         }
101         total_recv = rcount;
102     } else if (!(vrank % 2)) {
103       /* other non-leaf nodes, allocate temp buffer for data received from
104        * children, the most we need is half of the total data elements due
105        * to the property of binimoal tree */
106       tempbuf = (char*)smpi_get_tmp_sendbuffer(strue_extent + (scount * size - 1) * sextent);
107       if (NULL == tempbuf) {
108         err  = MPI_ERR_OTHER;
109         line = __LINE__;
110         goto err_hndl;
111       }
112
113       ptmp = tempbuf - slb;
114       /* local copy to tempbuf */
115       err = Datatype::copy(sbuf, scount, sdtype, ptmp, scount, sdtype);
116       if (MPI_SUCCESS != err) {
117         line = __LINE__;
118         goto err_hndl;
119       }
120
121       /* use sdtype,scount as rdtype,rdcount since they are ignored on
122        * non-root procs */
123       rdtype     = sdtype;
124       rcount     = scount;
125       rextent    = sextent;
126       total_recv = rcount;
127     } else {
128       /* leaf nodes, no temp buffer needed, use sdtype,scount as
129        * rdtype,rdcount since they are ignored on non-root procs */
130       ptmp       = (char*)sbuf;
131       total_recv = scount;
132     }
133
134     if (!(vrank % 2)) {
135       /* all non-leaf nodes recv from children */
136       for (i = 0; i < bmtree->tree_nextsize; i++) {
137         int mycount = 0, vkid;
138         /* figure out how much data I have to send to this child */
139         vkid    = (bmtree->tree_next[i] - root + size) % size;
140         mycount = vkid - vrank;
141         if (mycount > (size - vkid))
142           mycount = size - vkid;
143         mycount *= rcount;
144
145         XBT_DEBUG("smpi_coll_tuned_gather_ompi_binomial rank %d recv %d mycount = %d", rank, bmtree->tree_next[i],
146                   mycount);
147
148         Request::recv(ptmp + total_recv * rextent, mycount, rdtype, bmtree->tree_next[i], COLL_TAG_GATHER, comm,
149                       &status);
150
151         total_recv += mycount;
152       }
153     }
154
155     if (rank != root) {
156       /* all nodes except root send to parents */
157       XBT_DEBUG("smpi_coll_tuned_gather_ompi_binomial rank %d send %d count %d\n", rank, bmtree->tree_prev, total_recv);
158
159       Request::send(ptmp, total_recv, sdtype, bmtree->tree_prev, COLL_TAG_GATHER, comm);
160   }
161     if (rank == root) {
162       if (root != 0) {
163         /* rotate received data on root if root != 0 */
164         err = Datatype::copy(ptmp, rcount * (size - root), rdtype, (char*)rbuf + rextent * root * rcount,
165                              rcount * (size - root), rdtype);
166         if (MPI_SUCCESS != err) {
167           line = __LINE__;
168           goto err_hndl;
169         }
170
171         err = Datatype::copy(ptmp + rextent * rcount * (size - root), rcount * root, rdtype, (char*)rbuf, rcount * root,
172                              rdtype);
173         if (MPI_SUCCESS != err) {
174           line = __LINE__;
175           goto err_hndl;
176         }
177
178         smpi_free_tmp_buffer(tempbuf);
179       }
180     } else if (!(vrank % 2)) {
181       /* other non-leaf nodes */
182       smpi_free_tmp_buffer(tempbuf);
183     }
184     xbt_free(bmtree);
185     return MPI_SUCCESS;
186
187  err_hndl:
188     if (NULL != tempbuf)
189       smpi_free_tmp_buffer(tempbuf);
190
191     XBT_DEBUG("%s:%4d\tError occurred %d, rank %2d", __FILE__, line, err, rank);
192     return err;
193 }
194
195 /*
196  *  gather_intra_linear_sync
197  *
198  *  Function:  - synchronized gather operation with
199  *  Accepts:  - same arguments as MPI_Gather(), first segment size
200  *  Returns:  - MPI_SUCCESS or error code
201  */
202 int Coll_gather_ompi_linear_sync::gather(void *sbuf, int scount,
203                                          MPI_Datatype sdtype,
204                                          void *rbuf, int rcount,
205                                          MPI_Datatype rdtype,
206                                          int root,
207                                          MPI_Comm comm)
208 {
209     int i;
210     int ret, line;
211     int rank, size;
212     int first_segment_count;
213     size_t typelng;
214     MPI_Aint extent;
215     MPI_Aint lb;
216
217     int first_segment_size=0;
218     size = comm->size();
219     rank = comm->rank();
220
221     size_t dsize, block_size;
222     if (rank == root) {
223         dsize= rdtype->size();
224         block_size = dsize * rcount;
225     } else {
226         dsize=sdtype->size();
227         block_size = dsize * scount;
228     }
229
230      if (block_size > 92160){
231      first_segment_size = 32768;
232      }else{
233      first_segment_size = 1024;
234      }
235
236      XBT_DEBUG("smpi_coll_tuned_gather_ompi_linear_sync rank %d, segment %d", rank, first_segment_size);
237
238      if (rank != root) {
239        /* Non-root processes:
240           - receive zero byte message from the root,
241           - send the first segment of the data synchronously,
242           - send the second segment of the data.
243        */
244
245        typelng = sdtype->size();
246        sdtype->extent(&lb, &extent);
247        first_segment_count = scount;
248        COLL_TUNED_COMPUTED_SEGCOUNT((size_t)first_segment_size, typelng, first_segment_count);
249
250        Request::recv(sbuf, 0, MPI_BYTE, root, COLL_TAG_GATHER, comm, MPI_STATUS_IGNORE);
251
252        Request::send(sbuf, first_segment_count, sdtype, root, COLL_TAG_GATHER, comm);
253
254        Request::send((char*)sbuf + extent * first_segment_count, (scount - first_segment_count), sdtype, root,
255                      COLL_TAG_GATHER, comm);
256     }
257
258     else {
259       /* Root process,
260          - For every non-root node:
261    - post irecv for the first segment of the message
262    - send zero byte message to signal node to send the message
263    - post irecv for the second segment of the message
264    - wait for the first segment to complete
265          - Copy local data if necessary
266          - Waitall for all the second segments to complete.
267 */
268       char* ptmp;
269       MPI_Request *reqs = NULL, first_segment_req;
270       reqs              = (MPI_Request*)calloc(size, sizeof(MPI_Request));
271       if (NULL == reqs) {
272         ret  = -1;
273         line = __LINE__;
274         goto error_hndl; }
275
276         typelng=rdtype->size();
277         rdtype->extent(&lb, &extent);
278         first_segment_count = rcount;
279         COLL_TUNED_COMPUTED_SEGCOUNT( (size_t)first_segment_size, typelng,
280                                       first_segment_count );
281
282         for (i = 0; i < size; ++i) {
283             if (i == rank) {
284                 /* skip myself */
285                 reqs[i] = MPI_REQUEST_NULL;
286                 continue;
287             }
288
289             /* irecv for the first segment from i */
290             ptmp = (char*)rbuf + i * rcount * extent;
291             first_segment_req = Request::irecv(ptmp, first_segment_count, rdtype, i,
292                                      COLL_TAG_GATHER, comm
293                                      );
294
295             /* send sync message */
296             Request::send(rbuf, 0, MPI_BYTE, i,
297                                     COLL_TAG_GATHER,
298                                      comm);
299
300             /* irecv for the second segment */
301             ptmp = (char*)rbuf + (i * rcount + first_segment_count) * extent;
302             reqs[i]=Request::irecv(ptmp, (rcount - first_segment_count),
303                                      rdtype, i, COLL_TAG_GATHER, comm
304                                      );
305
306             /* wait on the first segment to complete */
307             Request::wait(&first_segment_req, MPI_STATUS_IGNORE);
308         }
309
310         /* copy local data if necessary */
311         if (MPI_IN_PLACE != sbuf) {
312             ret = Datatype::copy(sbuf, scount, sdtype,
313                                   (char*)rbuf + rank * rcount * extent,
314                                   rcount, rdtype);
315             if (ret != MPI_SUCCESS) { line = __LINE__; goto error_hndl; }
316         }
317
318         /* wait all second segments to complete */
319         ret = Request::waitall(size, reqs, MPI_STATUSES_IGNORE);
320         if (ret != MPI_SUCCESS) { line = __LINE__; goto error_hndl; }
321
322         free(reqs);
323     }
324
325     /* All done */
326
327     return MPI_SUCCESS;
328  error_hndl:
329     XBT_DEBUG(
330                    "ERROR_HNDL: node %d file %s line %d error %d\n",
331                    rank, __FILE__, line, ret );
332     return ret;
333 }
334
335 /*
336  * Linear functions are copied from the BASIC coll module
337  * they do not segment the message and are simple implementations
338  * but for some small number of nodes and/or small data sizes they
339  * are just as fast as tuned/tree based segmenting operations
340  * and as such may be selected by the decision functions
341  * These are copied into this module due to the way we select modules
342  * in V1. i.e. in V2 we will handle this differently and so will not
343  * have to duplicate code.
344  * JPG following the examples from other coll_tuned implementations. Dec06.
345  */
346
347 /* copied function (with appropriate renaming) starts here */
348 /*
349  *  gather_intra
350  *
351  *  Function:  - basic gather operation
352  *  Accepts:  - same arguments as MPI_Gather()
353  *  Returns:  - MPI_SUCCESS or error code
354  */
355 int Coll_gather_ompi_basic_linear::gather(void* sbuf, int scount, MPI_Datatype sdtype, void* rbuf, int rcount,
356                                           MPI_Datatype rdtype, int root, MPI_Comm comm)
357 {
358     int i;
359     int err;
360     int rank;
361     int size;
362     char *ptmp;
363     MPI_Aint incr;
364     MPI_Aint extent;
365     MPI_Aint lb;
366
367     size = comm->size();
368     rank = comm->rank();
369
370     /* Everyone but root sends data and returns. */
371     XBT_DEBUG("ompi_coll_tuned_gather_intra_basic_linear rank %d", rank);
372
373     if (rank != root) {
374         Request::send(sbuf, scount, sdtype, root,
375                                  COLL_TAG_GATHER,
376                                   comm);
377         return MPI_SUCCESS;
378     }
379
380     /* I am the root, loop receiving the data. */
381
382     rdtype->extent(&lb, &extent);
383     incr = extent * rcount;
384     for (i = 0, ptmp = (char *) rbuf; i < size; ++i, ptmp += incr) {
385         if (i == rank) {
386             if (MPI_IN_PLACE != sbuf) {
387                 err = Datatype::copy(sbuf, scount, sdtype,
388                                       ptmp, rcount, rdtype);
389             } else {
390                 err = MPI_SUCCESS;
391             }
392         } else {
393             Request::recv(ptmp, rcount, rdtype, i,
394                                     COLL_TAG_GATHER,
395                                     comm, MPI_STATUS_IGNORE);
396             err = MPI_SUCCESS;
397         }
398         if (MPI_SUCCESS != err) {
399             return err;
400         }
401     }
402
403     /* All done */
404
405     return MPI_SUCCESS;
406 }
407
408 }
409 }