Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' into actor-yield
[simgrid.git] / src / smpi / colls / allgather / allgather-NTSLR-NB.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.hpp"
8
9 namespace simgrid{
10 namespace smpi{
11
12
13 // Allgather-Non-Topoloty-Scecific-Logical-Ring algorithm
14 int
15 Coll_allgather_NTSLR_NB::allgather(void *sbuf, int scount, MPI_Datatype stype,
16                                    void *rbuf, int rcount, MPI_Datatype rtype,
17                                    MPI_Comm comm)
18 {
19   MPI_Aint rextent, sextent;
20   MPI_Status status, status2;
21   int i, to, from, rank, size;
22   int send_offset, recv_offset;
23   int tag = COLL_TAG_ALLGATHER;
24
25   rank = comm->rank();
26   size = comm->size();
27   rextent = rtype->get_extent();
28   sextent = stype->get_extent();
29   MPI_Request *rrequest_array;
30   MPI_Request *srequest_array;
31   rrequest_array = (MPI_Request *) xbt_malloc(size * sizeof(MPI_Request));
32   srequest_array = (MPI_Request *) xbt_malloc(size * sizeof(MPI_Request));
33
34   // irregular case use default MPI fucntions
35   if (scount * sextent != rcount * rextent) {
36     XBT_WARN("MPI_allgather_NTSLR_NB use default MPI_allgather.");
37     Coll_allgather_default::allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
38     return MPI_SUCCESS;
39   }
40
41   // topo non-specific
42   to = (rank + 1) % size;
43   from = (rank + size - 1) % size;
44
45   //copy a single segment from sbuf to rbuf
46   send_offset = rank * scount * sextent;
47
48   Request::sendrecv(sbuf, scount, stype, rank, tag,
49                (char *)rbuf + send_offset, rcount, rtype, rank, tag, comm, &status);
50
51
52   //start sending logical ring message
53   int increment = scount * sextent;
54
55   //post all irecv first
56   for (i = 0; i < size - 1; i++) {
57     recv_offset = ((rank - i - 1 + size) % size) * increment;
58     rrequest_array[i] = Request::irecv((char *)rbuf + recv_offset, rcount, rtype, from, tag + i, comm);
59   }
60
61
62   for (i = 0; i < size - 1; i++) {
63     send_offset = ((rank - i + size) % size) * increment;
64     srequest_array[i] = Request::isend((char *)rbuf + send_offset, scount, stype, to, tag + i, comm);
65     Request::wait(&rrequest_array[i], &status);
66     Request::wait(&srequest_array[i], &status2);
67   }
68
69   free(rrequest_array);
70   free(srequest_array);
71
72   return MPI_SUCCESS;
73 }
74
75 }
76 }