Logo AND Algorithmique Numérique Distribuée

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