Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
- MPI_Barrier() as a collective
[simgrid.git] / src / smpi / smpi_mpi.c
1
2
3 #include "private.h"
4 #include "smpi_coll_private.h"
5
6 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
7                                 "Logging specific to SMPI (mpi)");
8
9 int SMPI_MPI_Init(int *argc, char ***argv)
10 {
11   smpi_process_init(argc, argv);
12   smpi_bench_begin();
13   return MPI_SUCCESS;
14 }
15
16 int SMPI_MPI_Finalize()
17 {
18   smpi_bench_end();
19   smpi_process_finalize();
20   return MPI_SUCCESS;
21 }
22
23 // right now this just exits the current node, should send abort signal to all
24 // hosts in the communicator (TODO)
25 int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
26 {
27   smpi_exit(errorcode);
28   return 0;
29 }
30
31 int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
32 {
33   int retval = MPI_SUCCESS;
34
35   smpi_bench_end();
36
37   if (NULL == comm) {
38     retval = MPI_ERR_COMM;
39   } else if (NULL == size) {
40     retval = MPI_ERR_ARG;
41   } else {
42     *size = comm->size;
43   }
44
45   smpi_bench_begin();
46
47   return retval;
48 }
49
50 int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
51 {
52   int retval = MPI_SUCCESS;
53
54   smpi_bench_end();
55
56   if (NULL == comm) {
57     retval = MPI_ERR_COMM;
58   } else if (NULL == rank) {
59     retval = MPI_ERR_ARG;
60   } else {
61     *rank = smpi_mpi_comm_rank(comm);
62   }
63
64   smpi_bench_begin();
65
66   return retval;
67 }
68
69 int SMPI_MPI_Type_size(MPI_Datatype datatype, size_t * size)
70 {
71   int retval = MPI_SUCCESS;
72
73   smpi_bench_end();
74
75   if (NULL == datatype) {
76     retval = MPI_ERR_TYPE;
77   } else if (NULL == size) {
78     retval = MPI_ERR_ARG;
79   } else {
80     *size = datatype->size;
81   }
82
83   smpi_bench_begin();
84
85   return retval;
86 }
87
88 int SMPI_MPI_Barrier(MPI_Comm comm)
89 {
90   int retval = MPI_SUCCESS;
91   int arity=4;
92
93   smpi_bench_end();
94
95   if (NULL == comm) {
96     retval = MPI_ERR_COMM;
97   } else {
98
99     /*
100      * original implemantation:
101      * retval = smpi_mpi_barrier(comm);
102      * this one is unrealistic: it just cond_waits, means no time.
103      */
104      retval = nary_tree_barrier( comm, arity );
105   }
106
107   smpi_bench_begin();
108
109   return retval;
110 }
111
112 int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
113                    int tag, MPI_Comm comm, MPI_Request * request)
114 {
115   int retval = MPI_SUCCESS;
116
117   smpi_bench_end();
118
119   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
120                                request);
121   if (NULL != *request && MPI_SUCCESS == retval) {
122     retval = smpi_mpi_irecv(*request);
123   }
124
125   smpi_bench_begin();
126
127   return retval;
128 }
129
130 int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
131                   int tag, MPI_Comm comm, MPI_Status * status)
132 {
133   int retval = MPI_SUCCESS;
134   smpi_mpi_request_t request;
135
136   smpi_bench_end();
137
138   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
139                                &request);
140   if (NULL != request && MPI_SUCCESS == retval) {
141     retval = smpi_mpi_irecv(request);
142     if (MPI_SUCCESS == retval) {
143       retval = smpi_mpi_wait(request, status);
144     }
145     xbt_mallocator_release(smpi_global->request_mallocator, request);
146   }
147
148   smpi_bench_begin();
149
150   return retval;
151 }
152
153 int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
154                    int tag, MPI_Comm comm, MPI_Request * request)
155 {
156   int retval = MPI_SUCCESS;
157
158   smpi_bench_end();
159
160   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
161                                request);
162   if (NULL != *request && MPI_SUCCESS == retval) {
163     retval = smpi_mpi_isend(*request);
164   }
165
166   smpi_bench_begin();
167
168   return retval;
169 }
170
171 int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
172                   int tag, MPI_Comm comm)
173 {
174   int retval = MPI_SUCCESS;
175   smpi_mpi_request_t request;
176
177   smpi_bench_end();
178
179   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
180                                &request);
181   if (NULL != request && MPI_SUCCESS == retval) {
182     retval = smpi_mpi_isend(request);
183     if (MPI_SUCCESS == retval) {
184       smpi_mpi_wait(request, MPI_STATUS_IGNORE);
185     }
186     xbt_mallocator_release(smpi_global->request_mallocator, request);
187   }
188
189   smpi_bench_begin();
190
191   return retval;
192 }
193
194 /**
195  * MPI_Wait and friends
196  **/
197 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
198 {
199   return smpi_mpi_wait(*request, status);
200 }
201
202 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
203 {
204   return smpi_mpi_waitall(count, requests, status);
205 }
206
207 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
208                      MPI_Status status[])
209 {
210   return smpi_mpi_waitany(count, requests, index, status);
211 }
212
213 /**
214  * MPI_Bcast
215  **/
216
217 /**
218  * flat bcast 
219  **/
220 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
221 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
222                 MPI_Comm comm)
223 {
224         int rank;
225         int retval = MPI_SUCCESS;
226         smpi_mpi_request_t request;
227
228         rank = smpi_mpi_comm_rank(comm);
229         if (rank == root) {
230                 retval = smpi_create_request(buf, count, datatype, root,
231                                 (root + 1) % comm->size, 0, comm, &request);
232                 request->forward = comm->size - 1;
233                 smpi_mpi_isend(request);
234         } else {
235                 retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
236                                 0, comm, &request);
237                 smpi_mpi_irecv(request);
238         }
239
240         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
241         xbt_mallocator_release(smpi_global->request_mallocator, request);
242
243         return(retval);
244
245 }
246
247 /**
248  * Bcast user entry point
249  **/
250 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
251                    MPI_Comm comm)
252 {
253   int retval = MPI_SUCCESS;
254
255   smpi_bench_end();
256
257   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
258   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
259
260   smpi_bench_begin();
261
262   return retval;
263 }
264
265
266
267 #ifdef DEBUG_REDUCE
268 /**
269  * debugging helper function
270  **/
271 static void print_buffer_int(void *buf, int len, const char *msg, int rank)
272 {
273   int tmp, *v;
274   printf("**[%d] %s: ", rank, msg);
275   for (tmp = 0; tmp < len; tmp++) {
276     v = buf;
277     printf("[%d]", v[tmp]);
278   }
279   printf("\n");
280   free(msg);
281 }
282 #endif
283
284 /**
285  * MPI_Reduce
286  **/
287 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
288                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
289 {
290         int retval = MPI_SUCCESS;
291         int rank;
292         int size;
293         int i;
294         int tag = 0;
295         smpi_mpi_request_t *requests;
296         smpi_mpi_request_t request;
297
298         smpi_bench_end();
299
300         rank = smpi_mpi_comm_rank(comm);
301         size = comm->size;
302
303         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
304
305 #ifdef DEBUG_REDUCE
306                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
307 #endif
308                 retval =
309                         smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
310                                         &request);
311                 smpi_mpi_isend(request);
312                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
313                 xbt_mallocator_release(smpi_global->request_mallocator, request);
314
315         } else {
316                 // i am the ROOT: wait for all buffers by creating one request by sender
317                 int src;
318                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
319
320                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
321                 for (i = 0; i < size-1; i++) {
322                         // we need 1 buffer per request to store intermediate receptions
323                         tmpbufs[i] = xbt_malloc(count * datatype->size);
324                 }  
325                 // root: initiliaze recv buf with my own snd buf
326                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
327
328                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
329                 // since we should op values as soon as one receiving request matches.
330                 for (i = 0; i < size-1; i++) {
331                         // reminder: for smpi_create_request() the src is always the process sending.
332                         src = i < root ? i : i + 1;
333                         retval = smpi_create_request(tmpbufs[i], count, datatype,
334                                         src, root, tag, comm, &(requests[i]));
335                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
336                                 if (MPI_SUCCESS == retval) {
337                                         smpi_mpi_irecv(requests[i]);
338                                 }
339                         }
340                 }
341                 // now, wait for completion of all irecv's.
342                 for (i = 0; i < size-1; i++) {
343                         int index = MPI_UNDEFINED;
344                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
345 #ifdef DEBUG_REDUCE
346                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
347                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
348                                         rank);
349 #endif
350
351                         // arg 2 is modified
352                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
353 #ifdef DEBUG_REDUCE
354                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
355 #endif
356                         xbt_free(tmpbufs[index]);
357                         /* FIXME: with the following line, it  generates an
358                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
359                          */
360                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
361                 }
362                 xbt_free(requests);
363                 xbt_free(tmpbufs);
364         }
365         smpi_bench_begin();
366         return retval;
367 }
368
369 /**
370  * MPI_Allreduce
371  *
372  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
373  **/
374 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
375                          MPI_Op op, MPI_Comm comm );
376 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
377                          MPI_Op op, MPI_Comm comm )
378 {
379   int retval = MPI_SUCCESS;
380   int root=1;  // arbitrary choice
381
382   smpi_bench_end();
383
384   retval = SMPI_MPI_Reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
385   if (MPI_SUCCESS != retval)
386             return(retval);
387
388   retval = SMPI_MPI_Bcast( sendbuf, count, datatype, root, comm);
389   smpi_bench_begin();
390   return( retval );
391 }
392
393
394
395
396
397 // used by comm_split to sort ranks based on key values
398 int smpi_compare_rankkeys(const void *a, const void *b);
399 int smpi_compare_rankkeys(const void *a, const void *b)
400 {
401   int *x = (int *) a;
402   int *y = (int *) b;
403
404   if (x[1] < y[1])
405     return -1;
406
407   if (x[1] == y[1]) {
408     if (x[0] < y[0])
409       return -1;
410     if (x[0] == y[0])
411       return 0;
412     return 1;
413   }
414
415   return 1;
416 }
417
418 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
419                         MPI_Comm * comm_out)
420 {
421   int retval = MPI_SUCCESS;
422
423   int index, rank;
424   smpi_mpi_request_t request;
425   int colorkey[2];
426   smpi_mpi_status_t status;
427
428   smpi_bench_end();
429
430   // FIXME: need to test parameters
431
432   index = smpi_process_index();
433   rank = comm->index_to_rank_map[index];
434
435   // default output
436   comm_out = NULL;
437
438   // root node does most of the real work
439   if (0 == rank) {
440     int colormap[comm->size];
441     int keymap[comm->size];
442     int rankkeymap[comm->size * 2];
443     int i, j;
444     smpi_mpi_communicator_t tempcomm = NULL;
445     int count;
446     int indextmp;
447
448     colormap[0] = color;
449     keymap[0] = key;
450
451     // FIXME: use scatter/gather or similar instead of individual comms
452     for (i = 1; i < comm->size; i++) {
453       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
454                                    rank, MPI_ANY_TAG, comm, &request);
455       smpi_mpi_irecv(request);
456       smpi_mpi_wait(request, &status);
457       colormap[status.MPI_SOURCE] = colorkey[0];
458       keymap[status.MPI_SOURCE] = colorkey[1];
459       xbt_mallocator_release(smpi_global->request_mallocator, request);
460     }
461
462     for (i = 0; i < comm->size; i++) {
463       if (MPI_UNDEFINED == colormap[i]) {
464         continue;
465       }
466       // make a list of nodes with current color and sort by keys
467       count = 0;
468       for (j = i; j < comm->size; j++) {
469         if (colormap[i] == colormap[j]) {
470           colormap[j] = MPI_UNDEFINED;
471           rankkeymap[count * 2] = j;
472           rankkeymap[count * 2 + 1] = keymap[j];
473           count++;
474         }
475       }
476       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
477
478       // new communicator
479       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
480       tempcomm->barrier_count = 0;
481       tempcomm->size = count;
482       tempcomm->barrier_mutex = SIMIX_mutex_init();
483       tempcomm->barrier_cond = SIMIX_cond_init();
484       tempcomm->rank_to_index_map = xbt_new(int, count);
485       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
486       for (j = 0; j < smpi_global->process_count; j++) {
487         tempcomm->index_to_rank_map[j] = -1;
488       }
489       for (j = 0; j < count; j++) {
490         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
491         tempcomm->rank_to_index_map[j] = indextmp;
492         tempcomm->index_to_rank_map[indextmp] = j;
493       }
494       for (j = 0; j < count; j++) {
495         if (rankkeymap[j * 2]) {
496           retval = smpi_create_request(&j, 1, MPI_INT, 0,
497                                        rankkeymap[j * 2], 0, comm, &request);
498           request->data = tempcomm;
499           smpi_mpi_isend(request);
500           smpi_mpi_wait(request, &status);
501           xbt_mallocator_release(smpi_global->request_mallocator, request);
502         } else {
503           *comm_out = tempcomm;
504         }
505       }
506     }
507   } else {
508     colorkey[0] = color;
509     colorkey[1] = key;
510     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
511                                  &request);
512     smpi_mpi_isend(request);
513     smpi_mpi_wait(request, &status);
514     xbt_mallocator_release(smpi_global->request_mallocator, request);
515     if (MPI_UNDEFINED != color) {
516       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
517                                    &request);
518       smpi_mpi_irecv(request);
519       smpi_mpi_wait(request, &status);
520       *comm_out = request->data;
521     }
522   }
523
524   smpi_bench_begin();
525
526   return retval;
527 }
528
529 double SMPI_MPI_Wtime(void)
530 {
531   return (SIMIX_get_clock());
532 }
533
534