Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Use simgrid function instead of MPI in collectives
[simgrid.git] / src / smpi / colls / allgather-NTSLR-NB.c
1 #include "colls_private.h"
2
3 // Allgather-Non-Topoloty-Scecific-Logical-Ring algorithm
4 int
5 smpi_coll_tuned_allgather_NTSLR_NB(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, status2;
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   MPI_Request *rrequest_array;
20   MPI_Request *srequest_array;
21   rrequest_array = (MPI_Request *) xbt_malloc(size * sizeof(MPI_Request));
22   srequest_array = (MPI_Request *) xbt_malloc(size * sizeof(MPI_Request));
23
24   // irregular case use default MPI fucntions
25   if (scount * sextent != rcount * rextent) {
26     XBT_WARN("MPI_allgather_NTSLR_NB use default MPI_allgather.");  
27     smpi_mpi_allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
28     return MPI_SUCCESS;    
29   }
30
31   // topo non-specific
32   to = (rank + 1) % size;
33   from = (rank + size - 1) % size;
34
35   //copy a single segment from sbuf to rbuf
36   send_offset = rank * scount * sextent;
37
38   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
39                (char *)rbuf + send_offset, rcount, rtype, rank, tag, comm, &status);
40
41
42   //start sending logical ring message
43   int increment = scount * sextent;
44
45   //post all irecv first
46   for (i = 0; i < size - 1; i++) {
47     recv_offset = ((rank - i - 1 + size) % size) * increment;
48     rrequest_array[i] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount, rtype, from, tag + i, comm);
49   }
50
51
52   for (i = 0; i < size - 1; i++) {
53     send_offset = ((rank - i + size) % size) * increment;
54     srequest_array[i] = smpi_mpi_isend((char *)rbuf + send_offset, scount, stype, to, tag + i, comm);
55     smpi_mpi_wait(&rrequest_array[i], &status);
56     smpi_mpi_wait(&srequest_array[i], &status2);
57   }
58
59   free(rrequest_array);
60   free(srequest_array);
61
62   return MPI_SUCCESS;
63 }