Logo AND Algorithmique Numérique Distribuée

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