Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove another algorithm
[simgrid.git] / src / smpi / colls / allgather-lr.c
1 #include "colls.h"
2
3 // Allgather-Non-Topoloty-Scecific-Logical-Ring algorithm
4 int
5 smpi_coll_tuned_allgather_lr(void *sbuf, int scount, MPI_Datatype stype,
6                              void *rbuf, int rcount, MPI_Datatype rtype,
7                              MPI_Comm comm)
8 {
9   MPI_Aint rextent, sextent;
10   MPI_Status status;
11   int i, to, from, rank, size;
12   int send_offset, recv_offset;
13   int tag = 500;
14
15   MPI_Comm_rank(comm, &rank);
16   MPI_Comm_size(comm, &size);
17   MPI_Type_extent(rtype, &rextent);
18   MPI_Type_extent(stype, &sextent);
19
20   // irregular case use default MPI fucntions
21   if (scount * sextent != rcount * rextent)
22     MPI_Allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
23
24   // topo non-specific
25   to = (rank + 1) % size;
26   from = (rank + size - 1) % size;
27
28   //copy a single segment from sbuf to rbuf
29   send_offset = rank * scount * sextent;
30   MPI_Sendrecv(sbuf, scount, stype, rank, tag,
31                (char *) rbuf + send_offset, rcount, rtype, rank, tag,
32                comm, &status);
33
34   //start sending logical ring message
35   int increment = scount * sextent;
36   for (i = 0; i < size - 1; i++) {
37     send_offset = ((rank - i + size) % size) * increment;
38     recv_offset = ((rank - i - 1 + size) % size) * increment;
39     MPI_Sendrecv((char *) rbuf + send_offset, scount, stype, to, tag + i,
40                  (char *) rbuf + recv_offset, rcount, rtype, from, tag + i,
41                  comm, &status);
42   }
43
44   return MPI_SUCCESS;
45 }