Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
potential fixes
[simgrid.git] / src / smpi / colls / allreduce-smp-rsag.c
1 #include "colls_private.h"
2
3 /* change number of core per smp-node
4    we assume that number of core per process will be the same for all implementations */
5 #ifndef NUM_CORE
6 #define NUM_CORE 8
7 #endif
8
9 /*
10 This fucntion performs all-reduce operation as follow.
11 1) binomial_tree reduce inside each SMP node
12 2) reduce-scatter -inter between root of each SMP node
13 3) allgather - inter between root of each SMP node
14 4) binomial_tree bcast inside each SMP node
15 */
16 int smpi_coll_tuned_allreduce_smp_rsag(void *send_buf, void *recv_buf,
17                                        int count, MPI_Datatype dtype, MPI_Op op,
18                                        MPI_Comm comm)
19 {
20   int comm_size, rank;
21   void *tmp_buf;
22   int tag = COLL_TAG_ALLREDUCE;
23   int mask, src, dst;
24   MPI_Status status;
25   int num_core = NUM_CORE;
26   /*
27      #ifdef MPICH2_REDUCTION
28      MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
29      #else
30      MPI_User_function *uop;
31      struct MPIR_OP *op_ptr;
32      op_ptr = MPIR_ToPointer(op);
33      uop  = op_ptr->op;
34      #endif
35    */
36   comm_size = smpi_comm_size(comm);
37   rank = smpi_comm_rank(comm);
38   MPI_Aint extent;
39   extent = smpi_datatype_get_extent(dtype);
40   tmp_buf = (void *) xbt_malloc(count * extent);
41
42   int intra_rank, inter_rank;
43   intra_rank = rank % num_core;
44   inter_rank = rank / num_core;
45
46   //printf("node %d intra_rank = %d, inter_rank = %d\n", rank, intra_rank, inter_rank);
47
48   int inter_comm_size = (comm_size + num_core - 1) / num_core;
49
50   if (!rank) {
51     //printf("intra com size = %d\n",num_core);
52     //printf("inter com size = %d\n",inter_comm_size);
53   }
54
55
56   smpi_mpi_sendrecv(send_buf, count, dtype, rank, tag,
57                recv_buf, count, dtype, rank, tag, comm, &status);
58
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         smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
68         smpi_op_apply(op, tmp_buf, recv_buf, &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       smpi_mpi_send(recv_buf, 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
83   // INTER: reduce-scatter
84   if (intra_rank == 0) {
85     int send_offset, recv_offset;
86     int seg_count = count / inter_comm_size;
87     int to = ((inter_rank + 1) % inter_comm_size) * num_core;
88     int from =
89         ((inter_rank + inter_comm_size - 1) % inter_comm_size) * num_core;
90     int i;
91
92     //printf("node %d to %d from %d\n",rank,to,from);
93
94     for (i = 0; i < (inter_comm_size - 1); i++) {
95
96       send_offset =
97           ((inter_rank - 1 - i +
98             inter_comm_size) % inter_comm_size) * seg_count * extent;
99       recv_offset =
100           ((inter_rank - 2 - i +
101             inter_comm_size) % inter_comm_size) * seg_count * extent;
102
103       smpi_mpi_sendrecv((char *) recv_buf + send_offset, seg_count, dtype, to,
104                    tag + i, tmp_buf, seg_count, dtype, from, tag + i, comm,
105                    &status);
106
107       // result is in rbuf
108       smpi_op_apply(op, tmp_buf, (char *) recv_buf + recv_offset, &seg_count,
109                      &dtype);
110     }
111
112     // INTER: allgather
113     for (i = 0; i < (inter_comm_size - 1); i++) {
114
115       send_offset =
116           ((inter_rank - i +
117             inter_comm_size) % inter_comm_size) * seg_count * extent;
118       recv_offset =
119           ((inter_rank - 1 - i +
120             inter_comm_size) % inter_comm_size) * seg_count * extent;
121
122       smpi_mpi_sendrecv((char *) recv_buf + send_offset, seg_count, dtype, to,
123                    tag + i, (char *) recv_buf + recv_offset, seg_count, dtype,
124                    from, tag + i, comm, &status);
125
126     }
127   }
128
129
130
131   /*
132      // INTER_binomial_reduce
133
134      // only root node for each SMP
135      if (intra_rank == 0) {
136
137      mask = 1;
138      while (mask < inter_comm_size) {
139      if ((mask & inter_rank) == 0) {
140      src = (inter_rank | mask) * num_core;
141      if (src < comm_size) {
142      smpi_mpi_recv(tmp_buf, count, dtype, src, tag, comm, &status);
143      (* uop) (tmp_buf, recv_buf, &count, &dtype);
144      //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
145      }
146      }
147      else {
148      dst = (inter_rank & (~mask)) * num_core;
149      smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
150      //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
151      break;
152      }
153      mask <<=1;
154      }
155      }
156    */
157
158   /*
159      // INTER_binomial_bcast
160
161
162      if (intra_rank == 0) {
163      mask = 1;
164      while (mask < inter_comm_size) {
165      if (inter_rank & mask) {
166      src = (inter_rank - mask) * num_core;
167      //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
168      smpi_mpi_recv(recv_buf, count, dtype, src, tag, comm, &status);
169      break;
170      }
171      mask <<= 1;
172      }
173
174      mask >>= 1;
175      //printf("My rank = %d my mask = %d\n", rank,mask);
176
177      while (mask > 0) {
178      if (inter_rank < inter_comm_size) {
179      dst = (inter_rank + mask) * num_core;
180      if (dst < comm_size) {
181      //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
182      smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
183      }
184      }
185      mask >>= 1;
186      }
187      }
188    */
189
190
191   // INTRA_binomial_bcast
192
193   int num_core_in_current_smp = num_core;
194   if (inter_rank == (inter_comm_size - 1)) {
195     num_core_in_current_smp = comm_size - (inter_rank * num_core);
196   }
197   //  printf("Node %d num_core = %d\n",rank, num_core_in_current_smp);
198   mask = 1;
199   while (mask < num_core_in_current_smp) {
200     if (intra_rank & mask) {
201       src = (inter_rank * num_core) + (intra_rank - mask);
202       //printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
203       smpi_mpi_recv(recv_buf, count, dtype, src, tag, comm, &status);
204       break;
205     }
206     mask <<= 1;
207   }
208
209   mask >>= 1;
210   //printf("My rank = %d my mask = %d\n", rank,mask);
211
212   while (mask > 0) {
213     dst = (inter_rank * num_core) + (intra_rank + mask);
214     if (dst < comm_size) {
215       //printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
216       smpi_mpi_send(recv_buf, count, dtype, dst, tag, comm);
217     }
218     mask >>= 1;
219   }
220
221
222   free(tmp_buf);
223   return MPI_SUCCESS;
224 }