Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
0a3d12ae5ef27a22bce0389b4940691cf7171c01
[simgrid.git] / src / smpi / colls / allreduce-smp-rsag-rab.c
1 #include "colls_private.h"
2 /* 
3  * implemented by Pitch Patarasuk, 07/01/2007
4  */
5 //#include <star-reduction.c>
6
7 /* change number of core per smp-node
8    we assume that number of core per process will be the same for all implementations */
9 #ifndef NUM_CORE
10 #define NUM_CORE 8
11 #endif
12
13 /*
14 This fucntion performs all-reduce operation as follow.
15 1) binomial_tree reduce inside each SMP node
16 2) reduce-scatter -inter between root of each SMP node
17 3) allgather - inter between root of each SMP node
18 4) binomial_tree bcast inside each SMP node
19 */
20 int smpi_coll_tuned_allreduce_smp_rsag_rab(void *sbuf, void *rbuf, int count,
21                                            MPI_Datatype dtype, MPI_Op op,
22                                            MPI_Comm comm)
23 {
24   int comm_size, rank;
25   void *tmp_buf;
26   int tag = COLL_TAG_ALLREDUCE;
27   int mask, src, dst;
28   MPI_Status status;
29   int num_core = simcall_host_get_core(SIMIX_host_self());
30   // do we use the default one or the number of cores in the platform ?
31   // if the number of cores is one, the platform may be simulated with 1 node = 1 core
32   if (num_core == 1) num_core = NUM_CORE;
33
34   comm_size = smpi_comm_size(comm);
35
36   if((comm_size&(comm_size-1)))
37     THROWF(arg_error,0, "allreduce smp rsag rab algorithm can't be used with non power of two number of processes ! ");
38
39   rank = smpi_comm_rank(comm);
40   MPI_Aint extent;
41   extent = smpi_datatype_get_extent(dtype);
42   tmp_buf = (void *) xbt_malloc(count * extent);
43
44   int intra_rank, inter_rank;
45   intra_rank = rank % num_core;
46   inter_rank = rank / num_core;
47
48   //printf("node %d intra_rank = %d, inter_rank = %d\n", rank, intra_rank, inter_rank);
49
50   int inter_comm_size = (comm_size + num_core - 1) / num_core;
51
52   smpi_mpi_sendrecv(sbuf, count, dtype, rank, tag,
53                rbuf, count, dtype, rank, tag, comm, &status);
54
55   // SMP_binomial_reduce
56   mask = 1;
57   while (mask < num_core) {
58     if ((mask & intra_rank) == 0) {
59       src = (inter_rank * num_core) + (intra_rank | mask);
60       //      if (src < ((inter_rank + 1) * num_core)) {
61       if (src < comm_size) {
62         smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
63         smpi_op_apply(op, tmp_buf, rbuf, &count, &dtype);
64         //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
65       }
66     } else {
67
68       dst = (inter_rank * num_core) + (intra_rank & (~mask));
69       smpi_mpi_send(rbuf, count, dtype, dst, tag, comm);
70       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
71       break;
72     }
73     mask <<= 1;
74   }
75
76
77   // INTER: reduce-scatter
78   if (intra_rank == 0) {
79
80     int dst, base_offset, send_base_offset, recv_base_offset, recv_chunk;
81     int curr_count, i, recv_offset, send_offset;
82
83     // reduce-scatter
84
85     recv_chunk = extent * count / (comm_size / num_core);
86
87     mask = 1;
88     i = 0;
89     curr_count = count / 2;
90     int phase = 0;
91     base_offset = 0;
92     send_base_offset = 0;
93     recv_base_offset = 0;
94
95     while (mask < (comm_size / num_core)) {
96       dst = inter_rank ^ mask;
97
98       // compute offsets
99       send_base_offset = base_offset;
100
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       smpi_mpi_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       smpi_op_apply(op, 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       smpi_mpi_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       smpi_mpi_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       smpi_mpi_send(rbuf, count, dtype, dst, tag, comm);
205     }
206     mask >>= 1;
207   }
208
209
210   free(tmp_buf);
211   return MPI_SUCCESS;
212 }