Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
protect (hopefully) collective communication algorithms from abuse.
[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 = NUM_CORE;
30
31   comm_size = smpi_comm_size(comm);
32
33   if((comm_size&(comm_size-1)))
34     THROWF(arg_error,0, "allreduce smp rsag rab algorithm can't be used with non power of two number of processes ! ");
35
36   rank = smpi_comm_rank(comm);
37   MPI_Aint extent;
38   extent = smpi_datatype_get_extent(dtype);
39   tmp_buf = (void *) xbt_malloc(count * extent);
40
41   int intra_rank, inter_rank;
42   intra_rank = rank % num_core;
43   inter_rank = rank / num_core;
44
45   //printf("node %d intra_rank = %d, inter_rank = %d\n", rank, intra_rank, inter_rank);
46
47   int inter_comm_size = (comm_size + num_core - 1) / num_core;
48
49   smpi_mpi_sendrecv(sbuf, count, dtype, rank, tag,
50                rbuf, count, dtype, rank, tag, comm, &status);
51
52   // SMP_binomial_reduce
53   mask = 1;
54   while (mask < num_core) {
55     if ((mask & intra_rank) == 0) {
56       src = (inter_rank * num_core) + (intra_rank | mask);
57       //      if (src < ((inter_rank + 1) * num_core)) {
58       if (src < comm_size) {
59         smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
60         smpi_op_apply(op, tmp_buf, rbuf, &count, &dtype);
61         //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
62       }
63     } else {
64
65       dst = (inter_rank * num_core) + (intra_rank & (~mask));
66       smpi_mpi_send(rbuf, count, dtype, dst, tag, comm);
67       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
68       break;
69     }
70     mask <<= 1;
71   }
72
73
74   // INTER: reduce-scatter
75   if (intra_rank == 0) {
76
77     int dst, base_offset, send_base_offset, recv_base_offset, recv_chunk;
78     int curr_count, i, recv_offset, send_offset;
79
80     // reduce-scatter
81
82     recv_chunk = extent * count / (comm_size / num_core);
83
84     mask = 1;
85     i = 0;
86     curr_count = count / 2;
87     int phase = 0;
88     base_offset = 0;
89     send_base_offset = 0;
90     recv_base_offset = 0;
91
92     while (mask < (comm_size / num_core)) {
93       dst = inter_rank ^ mask;
94
95       // compute offsets
96       send_base_offset = base_offset;
97
98       // right-handside
99       if (inter_rank & mask) {
100         recv_base_offset = base_offset + curr_count;
101         send_base_offset = base_offset;
102         base_offset = recv_base_offset;
103       }
104       // left-handside
105       else {
106         recv_base_offset = base_offset;
107         send_base_offset = base_offset + curr_count;
108       }
109       send_offset = send_base_offset * extent;
110       recv_offset = recv_base_offset * extent;
111
112       //      if (rank==7)
113       //      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);
114
115       smpi_mpi_sendrecv((char *)rbuf + send_offset, curr_count, dtype, (dst * num_core), tag,
116                    tmp_buf, curr_count, dtype, (dst * num_core), tag,
117                    comm, &status);
118
119       smpi_op_apply(op, tmp_buf, (char *)rbuf + recv_offset, &curr_count, &dtype);
120
121       mask *= 2;
122       curr_count /= 2;
123       phase++;
124     }
125
126
127     // INTER: allgather
128
129     int size = (comm_size / num_core) / 2;
130     base_offset = 0;
131     mask = 1;
132     while (mask < (comm_size / num_core)) {
133       if (inter_rank & mask) {
134         base_offset += size;
135       }
136       mask <<= 1;
137       size /= 2;
138     }
139
140     curr_count *= 2;
141     mask >>= 1;
142     i = 1;
143     phase = 0;
144     while (mask >= 1) {
145       // destination pair for both send and recv
146       dst = inter_rank ^ mask;
147
148       // compute offsets
149       send_base_offset = base_offset;
150       if (inter_rank & mask) {
151         recv_base_offset = base_offset - i;
152         base_offset -= i;
153       } else {
154         recv_base_offset = base_offset + i;
155       }
156       send_offset = send_base_offset * recv_chunk;
157       recv_offset = recv_base_offset * recv_chunk;
158
159       //      if (rank==7)
160       //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);
161
162       smpi_mpi_sendrecv((char *)rbuf + send_offset, curr_count, dtype, (dst * num_core), tag,
163                    (char *)rbuf + recv_offset, curr_count, dtype, (dst * num_core), tag,
164                    comm, &status);
165
166
167       curr_count *= 2;
168       i *= 2;
169       mask >>= 1;
170       phase++;
171     }
172
173
174   }                             // INTER
175
176   // intra SMP binomial bcast
177
178   int num_core_in_current_smp = num_core;
179   if (inter_rank == (inter_comm_size - 1)) {
180     num_core_in_current_smp = comm_size - (inter_rank * num_core);
181   }
182   //  printf("Node %d num_core = %d\n",rank, num_core_in_current_smp);
183   mask = 1;
184   while (mask < num_core_in_current_smp) {
185     if (intra_rank & mask) {
186       src = (inter_rank * num_core) + (intra_rank - mask);
187       //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
188       smpi_mpi_recv(rbuf, count, dtype, src, tag, comm, &status);
189       break;
190     }
191     mask <<= 1;
192   }
193
194   mask >>= 1;
195   //printf("My rank = %d my mask = %d\n", rank,mask);
196
197   while (mask > 0) {
198     dst = (inter_rank * num_core) + (intra_rank + mask);
199     if (dst < comm_size) {
200       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
201       smpi_mpi_send(rbuf, count, dtype, dst, tag, comm);
202     }
203     mask >>= 1;
204   }
205
206
207   free(tmp_buf);
208   return MPI_SUCCESS;
209 }