Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
At least. ignore ignorable
[simgrid.git] / src / smpi / smpi_mpi.c
1 #include "private.h"
2
3 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi, "Logging specific to SMPI (mpi)");
4
5 int SMPI_MPI_Init(int *argc, char ***argv)
6 {
7         smpi_mpi_init();
8         smpi_bench_begin();
9         return MPI_SUCCESS;
10 }
11
12 int SMPI_MPI_Finalize()
13 {
14         smpi_bench_end();
15         smpi_mpi_finalize();
16         return MPI_SUCCESS;
17 }
18
19 // right now this just exits the current node, should send abort signal to all
20 // hosts in the communicator;
21 int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
22 {
23         smpi_exit(errorcode);
24         return 0;
25 }
26
27 int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
28 {
29         int retval = MPI_SUCCESS;
30
31         smpi_bench_end();
32
33         if (NULL == comm) {
34                 retval = MPI_ERR_COMM;
35         } else if (NULL == size) {
36                 retval = MPI_ERR_ARG;
37         } else {
38                 *size = comm->size;
39         }
40
41         smpi_bench_begin();
42
43         return retval;
44 }
45
46 int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
47 {
48         int retval = MPI_SUCCESS;
49
50         smpi_bench_end();
51
52         if (NULL == comm) {
53                 retval = MPI_ERR_COMM;
54         } else if (NULL == rank) {
55                 retval = MPI_ERR_ARG;
56         } else {
57                 *rank = smpi_mpi_comm_rank(comm);
58         }
59
60         smpi_bench_begin();
61
62         return retval;
63 }
64
65 int SMPI_MPI_Type_size(MPI_Datatype datatype, size_t *size)
66 {
67         int retval = MPI_SUCCESS;
68
69         smpi_bench_end();
70
71         if (NULL == datatype) {
72                 retval = MPI_ERR_TYPE;
73         } else if (NULL == size) {
74                 retval = MPI_ERR_ARG;
75         } else {
76                 *size = datatype->size;
77         }
78
79         smpi_bench_begin();
80
81         return retval;
82 }
83
84 int SMPI_MPI_Barrier(MPI_Comm comm)
85 {
86         int retval = MPI_SUCCESS;
87
88         smpi_bench_end();
89
90         if (NULL == comm) {
91                 retval = MPI_ERR_COMM;
92         } else {
93                 retval = smpi_mpi_barrier(comm);
94         }
95
96         smpi_bench_begin();
97
98         return retval;
99 }
100
101 int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
102     int tag, MPI_Comm comm, MPI_Request *request)
103 {
104         int retval = MPI_SUCCESS;
105
106         smpi_bench_end();
107
108         retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
109         request);
110         if (NULL != *request && MPI_SUCCESS == retval) {
111                 retval = smpi_mpi_irecv(*request);
112         }
113
114         smpi_bench_begin();
115
116         return retval;
117 }
118
119 int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
120     int tag, MPI_Comm comm, MPI_Status *status)
121 {
122         int retval = MPI_SUCCESS;
123         smpi_mpi_request_t request;
124
125         smpi_bench_end();
126
127         retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
128         &request);
129         if (NULL != request && MPI_SUCCESS == retval) {
130                 retval = smpi_mpi_irecv(request);
131                 if (MPI_SUCCESS == retval) {
132                         retval = smpi_mpi_wait(request, status);
133                 }
134                 xbt_mallocator_release(smpi_global->request_mallocator, request);
135         }
136
137         smpi_bench_begin();
138
139         return retval;
140 }
141
142 int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
143     int tag, MPI_Comm comm, MPI_Request *request)
144 {
145         int retval = MPI_SUCCESS;
146
147         smpi_bench_end();
148
149         retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
150         request);
151         if (NULL != *request && MPI_SUCCESS == retval) {
152                 retval = smpi_mpi_isend(*request);
153         }
154
155         smpi_bench_begin();
156
157         return retval;
158 }
159
160 int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
161     int tag, MPI_Comm comm)
162 {
163         int retval = MPI_SUCCESS;
164         smpi_mpi_request_t request;
165
166         smpi_bench_end();
167
168         retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
169         &request);
170         if (NULL != request && MPI_SUCCESS == retval) {
171                 retval = smpi_mpi_isend(request);
172                 if (MPI_SUCCESS == retval) {
173                         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
174                 }
175                 xbt_mallocator_release(smpi_global->request_mallocator, request);
176         }
177
178         smpi_bench_begin();
179
180         return retval;
181 }
182
183 int SMPI_MPI_Wait(MPI_Request *request, MPI_Status *status) {
184         return smpi_mpi_wait(*request, status);
185 }
186
187 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
188     MPI_Comm comm)
189 {
190
191         int retval = MPI_SUCCESS;
192         int rank;
193         smpi_mpi_request_t request;
194
195         smpi_bench_end();
196
197         rank = smpi_mpi_comm_rank(comm);
198
199         if (rank == root) {
200                 retval = smpi_create_request(buf, count, datatype, root,
201             (root + 1) % comm->size, 0, comm, &request);
202                 request->forward = comm->size - 1;
203                 smpi_mpi_isend(request);
204         } else {
205                 retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
206             0, comm, &request);
207                 smpi_mpi_irecv(request);
208         }
209
210         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
211         xbt_mallocator_release(smpi_global->request_mallocator, request);
212
213         smpi_bench_begin();
214
215         return retval;
216 }
217
218 // used by comm_split to sort ranks based on key values
219 int smpi_compare_rankkeys(const void *a, const void *b);
220 int smpi_compare_rankkeys(const void *a, const void *b) {
221     int *x = (int *)a;
222     int *y = (int *)b;
223
224     if (x[1] < y[1]) return -1;
225
226     if (x[1] == y[1]) {
227         if (x[0] < y[0]) return -1;
228         if (x[0] == y[0]) return 0;
229         return 1;
230     }
231
232     return 1;
233 }
234
235 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm *comm_out)
236 {
237         int retval = MPI_SUCCESS;
238
239         int index, rank;
240         smpi_mpi_request_t request;
241         int colorkey[2];
242         smpi_mpi_status_t status;
243
244         smpi_bench_end();
245
246         // FIXME: need to test parameters
247
248         index = smpi_host_index();
249         rank  = comm->index_to_rank_map[index];
250
251     // default output
252     comm_out = NULL;
253
254     // root node does most of the real work
255         if (0 == rank) {
256                 int colormap[comm->size];
257         int keymap[comm->size];
258         int rankkeymap[comm->size * 2];
259                 int i, j;
260                 smpi_mpi_communicator_t tempcomm = NULL;
261                 int count;
262                 int indextmp;
263
264                 colormap[0] = color;
265         keymap[0]   = key;
266
267                 // FIXME: use scatter/gather or similar instead of individual comms
268                 for (i = 1; i < comm->size; i++) {
269                         retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
270                     rank, MPI_ANY_TAG, comm, &request);
271                         smpi_mpi_irecv(request);
272                         smpi_mpi_wait(request, &status);
273                         colormap[status.MPI_SOURCE] = colorkey[0];
274                         keymap[status.MPI_SOURCE]   = colorkey[1];
275                         xbt_mallocator_release(smpi_global->request_mallocator, request);
276                 }
277
278                 for (i = 0; i < comm->size; i++) {
279                         if (MPI_UNDEFINED == colormap[i]) {
280                                 continue;
281                         }
282
283             // make a list of nodes with current color and sort by keys
284                         count = 0;
285                         for (j = i; j < comm->size; j++) {
286                                 if(colormap[i] == colormap[j]) {
287                                         colormap[j] = MPI_UNDEFINED;
288                                         rankkeymap[count * 2]     = j;
289                                         rankkeymap[count * 2 + 1] = keymap[j];
290                                         count++;
291                                 }
292                         }
293             qsort(rankkeymap, count, sizeof(int) * 2,
294                 &smpi_compare_rankkeys);
295
296             // new communicator
297                         tempcomm                    = xbt_new(s_smpi_mpi_communicator_t, 1);
298                         tempcomm->barrier_count     = 0;
299                         tempcomm->size              = count;
300                         tempcomm->barrier_mutex     = SIMIX_mutex_init();
301                         tempcomm->barrier_cond      = SIMIX_cond_init();
302                         tempcomm->rank_to_index_map = xbt_new(int, count);
303                         tempcomm->index_to_rank_map = xbt_new(int, smpi_global->host_count);
304                         for (j = 0; j < smpi_global->host_count; j++) {
305                                 tempcomm->index_to_rank_map[j] = -1;
306                         }
307                         for (j = 0; j < count; j++) {
308                                 indextmp = comm->rank_to_index_map[rankkeymap[j*2]];
309                                 tempcomm->rank_to_index_map[j]        = indextmp;
310                                 tempcomm->index_to_rank_map[indextmp] = j;
311                         }
312                         for (j = 0; j < count; j++) {
313                                 if (rankkeymap[j*2]) {
314                                         retval = smpi_create_request(&j, 1, MPI_INT, 0,
315                         rankkeymap[j*2], 0, comm, &request);
316                                         request->data = tempcomm;
317                                         smpi_mpi_isend(request);
318                                         smpi_mpi_wait(request, &status);
319                                         xbt_mallocator_release(smpi_global->request_mallocator,
320                         request);
321                                 } else {
322                                         *comm_out = tempcomm;
323                                 }
324                         }
325                 }
326         } else {
327                 colorkey[0] = color;
328                 colorkey[1] = key;
329                 retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
330             &request);
331                 smpi_mpi_isend(request);
332                 smpi_mpi_wait(request, &status);
333                 xbt_mallocator_release(smpi_global->request_mallocator, request);
334         if (MPI_UNDEFINED != color) {
335                     retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
336                 &request);
337                     smpi_mpi_irecv(request);
338                     smpi_mpi_wait(request, &status);
339                     *comm_out = request->data;
340         }
341         }
342
343         smpi_bench_begin();
344
345         return retval;
346 }