Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add scatter algos from ompi
[simgrid.git] / src / smpi / colls / allgather-NTSLR.c
1 #include "colls_private.h"
2
3 // Allgather-Non-Topoloty-Scecific-Logical-Ring algorithm
4 int
5 smpi_coll_tuned_allgather_NTSLR(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   rank = smpi_comm_rank(comm);
16   size = smpi_comm_size(comm);
17   rextent = smpi_datatype_get_extent(rtype);
18   sextent = smpi_datatype_get_extent(stype);
19
20   // irregular case use default MPI fucntions
21   if (scount * sextent != rcount * rextent) {
22     XBT_WARN("MPI_allgather_NTSLR use default MPI_allgather.");  
23     smpi_mpi_allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
24     return MPI_SUCCESS;    
25   }
26
27   // topo non-specific
28   to = (rank + 1) % size;
29   from = (rank + size - 1) % size;
30
31   //copy a single segment from sbuf to rbuf
32   send_offset = rank * scount * sextent;
33
34   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
35                (char *)rbuf + send_offset, rcount, rtype, rank, tag,
36                comm, &status);
37
38
39   //start sending logical ring message
40   int increment = scount * sextent;
41   for (i = 0; i < size - 1; i++) {
42     send_offset = ((rank - i + size) % size) * increment;
43     recv_offset = ((rank - i - 1 + size) % size) * increment;
44     smpi_mpi_sendrecv((char *) rbuf + send_offset, scount, stype, to, tag + i,
45                  (char *) rbuf + recv_offset, rcount, rtype, from, tag + i,
46                  comm, &status);
47   }
48
49   return MPI_SUCCESS;
50 }