Logo AND Algorithmique Numérique Distribuée

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