Logo AND Algorithmique Numérique Distribuée

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