Logo AND Algorithmique Numérique Distribuée

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