Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add new entry in Release_Notes.
[simgrid.git] / src / smpi / colls / allreduce / allreduce-rdb.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 //#include <star-reduction.c>
9 namespace simgrid::smpi {
10 int allreduce__rdb(const void *sbuff, void *rbuff, int count,
11                    MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)
12 {
13   int nprocs, rank, tag = COLL_TAG_ALLREDUCE;
14   int mask, dst, pof2, newrank, rem, newdst;
15   MPI_Aint extent, lb;
16   MPI_Status status;
17   /*
18      #ifdef MPICH2_REDUCTION
19      MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
20      #else
21      MPI_User_function *uop;
22      MPIR_OP *op_ptr;
23      op_ptr = MPIR_ToPointer(op);
24      uop  = op_ptr->op;
25      #endif
26    */
27   nprocs=comm->size();
28   rank=comm->rank();
29
30   dtype->extent(&lb, &extent);
31   unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
32
33   Request::sendrecv(sbuff, count, dtype, rank, 500,
34                rbuff, count, dtype, rank, 500, comm, &status);
35
36   // find nearest power-of-two less than or equal to comm_size
37   pof2 = 1;
38   while (pof2 <= nprocs)
39     pof2 <<= 1;
40   pof2 >>= 1;
41
42   rem = nprocs - pof2;
43
44   // In the non-power-of-two case, all even-numbered
45   // processes of rank < 2*rem send their data to
46   // (rank+1). These even-numbered processes no longer
47   // participate in the algorithm until the very end. The
48   // remaining processes form a nice power-of-two.
49
50   if (rank < 2 * rem) {
51     // even
52     if (rank % 2 == 0) {
53
54       Request::send(rbuff, count, dtype, rank + 1, tag, comm);
55
56       // temporarily set the rank to -1 so that this
57       // process does not participate in recursive
58       // doubling
59       newrank = -1;
60     } else                      // odd
61     {
62       Request::recv(tmp_buf, count, dtype, rank - 1, tag, comm, &status);
63       // do the reduction on received data. since the
64       // ordering is right, it doesn't matter whether
65       // the operation is commutative or not.
66       if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuff, &count, dtype);
67
68       // change the rank
69       newrank = rank / 2;
70     }
71   }
72
73   else                          // rank >= 2 * rem
74     newrank = rank - rem;
75
76   // If op is user-defined or count is less than pof2, use
77   // recursive doubling algorithm. Otherwise do a reduce-scatter
78   // followed by allgather. (If op is user-defined,
79   // derived datatypes are allowed and the user could pass basic
80   // datatypes on one process and derived on another as long as
81   // the type maps are the same. Breaking up derived
82   // datatypes to do the reduce-scatter is tricky, therefore
83   // using recursive doubling in that case.)
84
85   if (newrank != -1) {
86     mask = 0x1;
87     while (mask < pof2) {
88       newdst = newrank ^ mask;
89       // find real rank of dest
90       dst = (newdst < rem) ? newdst * 2 + 1 : newdst + rem;
91
92       // Send the most current data, which is in recvbuf. Recv
93       // into tmp_buf
94       Request::sendrecv(rbuff, count, dtype, dst, tag, tmp_buf, count, dtype,
95                    dst, tag, comm, &status);
96
97       // tmp_buf contains data received in this step.
98       // recvbuf contains data accumulated so far
99
100       // op is commutative OR the order is already right
101       // we assume it is commutative op
102       //      if (op -> op_commute  || (dst < rank))
103       if ((dst < rank)) {
104         if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuff, &count, dtype);
105       } else                    // op is noncommutative and the order is not right
106       {
107         if(op!=MPI_OP_NULL) op->apply( rbuff, tmp_buf, &count, dtype);
108
109         // copy result back into recvbuf
110         Request::sendrecv(tmp_buf, count, dtype, rank, tag, rbuff, count,
111                      dtype, rank, tag, comm, &status);
112       }
113       mask <<= 1;
114     }
115   }
116   // In the non-power-of-two case, all odd-numbered processes of
117   // rank < 2 * rem send the result to (rank-1), the ranks who didn't
118   // participate above.
119
120   if (rank < 2 * rem) {
121     if (rank % 2)               // odd
122       Request::send(rbuff, count, dtype, rank - 1, tag, comm);
123     else                        // even
124       Request::recv(rbuff, count, dtype, rank + 1, tag, comm, &status);
125   }
126
127   smpi_free_tmp_buffer(tmp_buf);
128   return MPI_SUCCESS;
129 }
130 } // namespace simgrid::smpi