Logo AND Algorithmique Numérique Distribuée

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