Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
MPI_Scatter() ok
[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, 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 static void print_buffer_double(void *buf, int len, char *msg, int rank)
283 {
284   int tmp;
285   double *v;
286   printf("**[%d] %s: ", rank, msg);
287   for (tmp = 0; tmp < len; tmp++) {
288     v = buf;
289     printf("[%lf]", v[tmp]);
290   }
291   printf("\n");
292   free(msg);
293 }
294
295
296 //#endif
297 /**
298  * MPI_Reduce
299  **/
300 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
301                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
302 {
303         int retval = MPI_SUCCESS;
304         int rank;
305         int size;
306         int i;
307         int tag = 0;
308         smpi_mpi_request_t *requests;
309         smpi_mpi_request_t request;
310
311         smpi_bench_end();
312
313         rank = smpi_mpi_comm_rank(comm);
314         size = comm->size;
315
316         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
317
318 #ifdef DEBUG_REDUCE
319                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
320 #endif
321                 retval =
322                         smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
323                                         &request);
324                 smpi_mpi_isend(request);
325                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
326                 xbt_mallocator_release(smpi_global->request_mallocator, request);
327
328         } else {
329                 // i am the ROOT: wait for all buffers by creating one request by sender
330                 int src;
331                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
332
333                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
334                 for (i = 0; i < size-1; i++) {
335                         // we need 1 buffer per request to store intermediate receptions
336                         tmpbufs[i] = xbt_malloc(count * datatype->size);
337                 }  
338                 // root: initiliaze recv buf with my own snd buf
339                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
340
341                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
342                 // since we should op values as soon as one receiving request matches.
343                 for (i = 0; i < size-1; i++) {
344                         // reminder: for smpi_create_request() the src is always the process sending.
345                         src = i < root ? i : i + 1;
346                         retval = smpi_create_request(tmpbufs[i], count, datatype,
347                                         src, root, tag, comm, &(requests[i]));
348                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
349                                 if (MPI_SUCCESS == retval) {
350                                         smpi_mpi_irecv(requests[i]);
351                                 }
352                         }
353                 }
354                 // now, wait for completion of all irecv's.
355                 for (i = 0; i < size-1; i++) {
356                         int index = MPI_UNDEFINED;
357                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
358 #ifdef DEBUG_REDUCE
359                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
360                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
361                                         rank);
362 #endif
363
364                         // arg 2 is modified
365                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
366 #ifdef DEBUG_REDUCE
367                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
368 #endif
369                         xbt_free(tmpbufs[index]);
370                         /* FIXME: with the following line, it  generates an
371                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
372                          */
373                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
374                 }
375                 xbt_free(requests);
376                 xbt_free(tmpbufs);
377         }
378         smpi_bench_begin();
379         return retval;
380 }
381
382 /**
383  * MPI_Allreduce
384  *
385  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
386  **/
387 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
388                          MPI_Op op, MPI_Comm comm );
389 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
390                          MPI_Op op, MPI_Comm comm )
391 {
392   int retval = MPI_SUCCESS;
393   int root=1;  // arbitrary choice
394
395   smpi_bench_end();
396
397   retval = SMPI_MPI_Reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
398   if (MPI_SUCCESS != retval)
399             return(retval);
400
401   retval = SMPI_MPI_Bcast( sendbuf, count, datatype, root, comm);
402   smpi_bench_begin();
403   return( retval );
404 }
405
406
407 /**
408  * MPI_Scatter user entry point
409  **/
410 //int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
411 //                       void *recvbuf, int recvcount, MPI_Datatype recvtype,int root,
412 //                      MPI_Comm comm);
413 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
414                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
415                            int root, MPI_Comm comm)
416 {
417   int retval = MPI_SUCCESS;
418   int i;
419   int cnt=0;  
420   int rank;
421   int tag=0;
422   char *cptr;  // to manipulate the void * buffers
423   smpi_mpi_request_t *requests;
424   smpi_mpi_request_t request;
425   smpi_mpi_status_t status;
426
427
428   smpi_bench_end();
429
430   rank = smpi_mpi_comm_rank(comm);
431
432   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
433   if (rank == root) {
434           // i am the root: distribute my sendbuf
435           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
436           cptr = sendbuf;
437           for (i=0; i < comm->size; i++) {
438                   if ( i!=root ) { // send to processes ...
439
440                           retval = smpi_create_request((void *)cptr, sendcount, 
441                                           datatype, root, i, tag, comm, &(requests[cnt]));
442                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
443                                   if (MPI_SUCCESS == retval) {
444                                           smpi_mpi_isend(requests[cnt]);
445                                   }
446                                   }
447                                   cnt++;
448                         } 
449                         else { // ... except if it's me.
450                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
451                         }
452                   cptr += sendcount*datatype->size;
453             }
454             for(i=0; i<cnt; i++) { // wait for send to complete
455                             /* FIXME: waitall() should be slightly better */
456                             smpi_mpi_wait(requests[i], &status);
457                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
458
459             }
460   } 
461   else {  // i am a non-root process: wait data from the root
462             retval = smpi_create_request(recvbuf,recvcount, 
463                                   recvtype, root, rank, tag, comm, &request);
464             if (NULL != request && MPI_SUCCESS == retval) {
465                         if (MPI_SUCCESS == retval) {
466                                   smpi_mpi_irecv(request);
467                         }
468             }
469             smpi_mpi_wait(request, &status);
470             xbt_mallocator_release(smpi_global->request_mallocator, request);
471   }
472   xbt_free(requests);
473
474   smpi_bench_begin();
475
476   return retval;
477 }
478
479
480
481
482
483
484 // used by comm_split to sort ranks based on key values
485 int smpi_compare_rankkeys(const void *a, const void *b);
486 int smpi_compare_rankkeys(const void *a, const void *b)
487 {
488   int *x = (int *) a;
489   int *y = (int *) b;
490
491   if (x[1] < y[1])
492     return -1;
493
494   if (x[1] == y[1]) {
495     if (x[0] < y[0])
496       return -1;
497     if (x[0] == y[0])
498       return 0;
499     return 1;
500   }
501
502   return 1;
503 }
504
505 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
506                         MPI_Comm * comm_out)
507 {
508   int retval = MPI_SUCCESS;
509
510   int index, rank;
511   smpi_mpi_request_t request;
512   int colorkey[2];
513   smpi_mpi_status_t status;
514
515   smpi_bench_end();
516
517   // FIXME: need to test parameters
518
519   index = smpi_process_index();
520   rank = comm->index_to_rank_map[index];
521
522   // default output
523   comm_out = NULL;
524
525   // root node does most of the real work
526   if (0 == rank) {
527     int colormap[comm->size];
528     int keymap[comm->size];
529     int rankkeymap[comm->size * 2];
530     int i, j;
531     smpi_mpi_communicator_t tempcomm = NULL;
532     int count;
533     int indextmp;
534
535     colormap[0] = color;
536     keymap[0] = key;
537
538     // FIXME: use scatter/gather or similar instead of individual comms
539     for (i = 1; i < comm->size; i++) {
540       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
541                                    rank, MPI_ANY_TAG, comm, &request);
542       smpi_mpi_irecv(request);
543       smpi_mpi_wait(request, &status);
544       colormap[status.MPI_SOURCE] = colorkey[0];
545       keymap[status.MPI_SOURCE] = colorkey[1];
546       xbt_mallocator_release(smpi_global->request_mallocator, request);
547     }
548
549     for (i = 0; i < comm->size; i++) {
550       if (MPI_UNDEFINED == colormap[i]) {
551         continue;
552       }
553       // make a list of nodes with current color and sort by keys
554       count = 0;
555       for (j = i; j < comm->size; j++) {
556         if (colormap[i] == colormap[j]) {
557           colormap[j] = MPI_UNDEFINED;
558           rankkeymap[count * 2] = j;
559           rankkeymap[count * 2 + 1] = keymap[j];
560           count++;
561         }
562       }
563       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
564
565       // new communicator
566       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
567       tempcomm->barrier_count = 0;
568       tempcomm->size = count;
569       tempcomm->barrier_mutex = SIMIX_mutex_init();
570       tempcomm->barrier_cond = SIMIX_cond_init();
571       tempcomm->rank_to_index_map = xbt_new(int, count);
572       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
573       for (j = 0; j < smpi_global->process_count; j++) {
574         tempcomm->index_to_rank_map[j] = -1;
575       }
576       for (j = 0; j < count; j++) {
577         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
578         tempcomm->rank_to_index_map[j] = indextmp;
579         tempcomm->index_to_rank_map[indextmp] = j;
580       }
581       for (j = 0; j < count; j++) {
582         if (rankkeymap[j * 2]) {
583           retval = smpi_create_request(&j, 1, MPI_INT, 0,
584                                        rankkeymap[j * 2], 0, comm, &request);
585           request->data = tempcomm;
586           smpi_mpi_isend(request);
587           smpi_mpi_wait(request, &status);
588           xbt_mallocator_release(smpi_global->request_mallocator, request);
589         } else {
590           *comm_out = tempcomm;
591         }
592       }
593     }
594   } else {
595     colorkey[0] = color;
596     colorkey[1] = key;
597     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
598                                  &request);
599     smpi_mpi_isend(request);
600     smpi_mpi_wait(request, &status);
601     xbt_mallocator_release(smpi_global->request_mallocator, request);
602     if (MPI_UNDEFINED != color) {
603       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
604                                    &request);
605       smpi_mpi_irecv(request);
606       smpi_mpi_wait(request, &status);
607       *comm_out = request->data;
608     }
609   }
610
611   smpi_bench_begin();
612
613   return retval;
614 }
615
616 double SMPI_MPI_Wtime(void)
617 {
618   return (SIMIX_get_clock());
619 }
620
621