Logo AND Algorithmique Numérique Distribuée

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