Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
ad9c9bdfc90f8493e4331584b03acfdbc8a670b3
[simgrid.git] / src / smpi / colls / allreduce / allreduce-smp-rsag-rab.cpp
1 /* Copyright (c) 2013-2019. 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 /*
8  * implemented by Pitch Patarasuk, 07/01/2007
9  */
10 #include "../colls_private.hpp"
11 //#include <star-reduction.c>
12
13
14 /*
15 This fucntion performs all-reduce operation as follow.
16 1) binomial_tree reduce inside each SMP node
17 2) reduce-scatter -inter between root of each SMP node
18 3) allgather - inter between root of each SMP node
19 4) binomial_tree bcast inside each SMP node
20 */
21 namespace simgrid{
22 namespace smpi{
23 int Coll_allreduce_smp_rsag_rab::allreduce(const void *sbuf, void *rbuf, int count,
24                                            MPI_Datatype dtype, MPI_Op op,
25                                            MPI_Comm comm)
26 {
27   int comm_size, rank;
28   int tag = COLL_TAG_ALLREDUCE;
29   int mask, src, dst;
30   MPI_Status status;
31   if(comm->get_leaders_comm()==MPI_COMM_NULL){
32     comm->init_smp();
33   }
34   int num_core=1;
35   if (comm->is_uniform()){
36     num_core = comm->get_intra_comm()->size();
37   }
38
39   comm_size = comm->size();
40
41   if((comm_size&(comm_size-1)))
42     THROWF(arg_error,0, "allreduce smp rsag rab algorithm can't be used with non power of two number of processes ! ");
43
44   rank = comm->rank();
45   MPI_Aint extent;
46   extent = dtype->get_extent();
47   unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
48
49   int intra_rank, inter_rank;
50   intra_rank = rank % num_core;
51   inter_rank = rank / num_core;
52
53   //printf("node %d intra_rank = %d, inter_rank = %d\n", rank, intra_rank, inter_rank);
54
55   int inter_comm_size = (comm_size + num_core - 1) / num_core;
56
57   Request::sendrecv(sbuf, count, dtype, rank, tag,
58                rbuf, count, dtype, rank, tag, comm, &status);
59
60   // SMP_binomial_reduce
61   mask = 1;
62   while (mask < num_core) {
63     if ((mask & intra_rank) == 0) {
64       src = (inter_rank * num_core) + (intra_rank | mask);
65       //      if (src < ((inter_rank + 1) * num_core)) {
66       if (src < comm_size) {
67         Request::recv(tmp_buf, count, dtype, src, tag, comm, &status);
68         if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuf, &count, dtype);
69         //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
70       }
71     } else {
72
73       dst = (inter_rank * num_core) + (intra_rank & (~mask));
74       Request::send(rbuf, count, dtype, dst, tag, comm);
75       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
76       break;
77     }
78     mask <<= 1;
79   }
80
81
82   // INTER: reduce-scatter
83   if (intra_rank == 0) {
84
85     int dst, base_offset, send_base_offset, recv_base_offset, recv_chunk;
86     int curr_count, i, recv_offset, send_offset;
87
88     // reduce-scatter
89
90     recv_chunk = extent * count / (comm_size / num_core);
91
92     mask = 1;
93     curr_count = count / 2;
94     int phase = 0;
95     base_offset = 0;
96
97     while (mask < (comm_size / num_core)) {
98       dst = inter_rank ^ mask;
99
100       // compute offsets
101       // right-handside
102       if (inter_rank & mask) {
103         recv_base_offset = base_offset + curr_count;
104         send_base_offset = base_offset;
105         base_offset = recv_base_offset;
106       }
107       // left-handside
108       else {
109         recv_base_offset = base_offset;
110         send_base_offset = base_offset + curr_count;
111       }
112       send_offset = send_base_offset * extent;
113       recv_offset = recv_base_offset * extent;
114
115       //      if (rank==7)
116       //      printf("node %d send to %d in phase %d s_offset = %d r_offset = %d count = %d\n",rank,dst,phase, send_offset, recv_offset, curr_count);
117
118       Request::sendrecv((char *)rbuf + send_offset, curr_count, dtype, (dst * num_core), tag,
119                    tmp_buf, curr_count, dtype, (dst * num_core), tag,
120                    comm, &status);
121
122       if(op!=MPI_OP_NULL) op->apply( tmp_buf, (char *)rbuf + recv_offset, &curr_count, dtype);
123
124       mask *= 2;
125       curr_count /= 2;
126       phase++;
127     }
128
129
130     // INTER: allgather
131
132     int size = (comm_size / num_core) / 2;
133     base_offset = 0;
134     mask = 1;
135     while (mask < (comm_size / num_core)) {
136       if (inter_rank & mask) {
137         base_offset += size;
138       }
139       mask <<= 1;
140       size /= 2;
141     }
142
143     curr_count *= 2;
144     mask >>= 1;
145     i = 1;
146     phase = 0;
147     while (mask >= 1) {
148       // destination pair for both send and recv
149       dst = inter_rank ^ mask;
150
151       // compute offsets
152       send_base_offset = base_offset;
153       if (inter_rank & mask) {
154         recv_base_offset = base_offset - i;
155         base_offset -= i;
156       } else {
157         recv_base_offset = base_offset + i;
158       }
159       send_offset = send_base_offset * recv_chunk;
160       recv_offset = recv_base_offset * recv_chunk;
161
162       //      if (rank==7)
163       //printf("node %d send to %d in phase %d s_offset = %d r_offset = %d count = %d\n",rank,dst,phase, send_offset, recv_offset, curr_count);
164
165       Request::sendrecv((char *)rbuf + send_offset, curr_count, dtype, (dst * num_core), tag,
166                    (char *)rbuf + recv_offset, curr_count, dtype, (dst * num_core), tag,
167                    comm, &status);
168
169
170       curr_count *= 2;
171       i *= 2;
172       mask >>= 1;
173       phase++;
174     }
175
176
177   }                             // INTER
178
179   // intra SMP binomial bcast
180
181   int num_core_in_current_smp = num_core;
182   if (inter_rank == (inter_comm_size - 1)) {
183     num_core_in_current_smp = comm_size - (inter_rank * num_core);
184   }
185   //  printf("Node %d num_core = %d\n",rank, num_core_in_current_smp);
186   mask = 1;
187   while (mask < num_core_in_current_smp) {
188     if (intra_rank & mask) {
189       src = (inter_rank * num_core) + (intra_rank - mask);
190       //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
191       Request::recv(rbuf, count, dtype, src, tag, comm, &status);
192       break;
193     }
194     mask <<= 1;
195   }
196
197   mask >>= 1;
198   //printf("My rank = %d my mask = %d\n", rank,mask);
199
200   while (mask > 0) {
201     dst = (inter_rank * num_core) + (intra_rank + mask);
202     if (dst < comm_size) {
203       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
204       Request::send(rbuf, count, dtype, dst, tag, comm);
205     }
206     mask >>= 1;
207   }
208
209
210   smpi_free_tmp_buffer(tmp_buf);
211   return MPI_SUCCESS;
212 }
213 }
214 }