Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
9cf5d01ddc5c07c3527394120e74c3759660eb2d
[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_Sendrecv
196  **/
197 int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
198                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
199                     MPI_Comm comm, MPI_Status *status)
200 {
201 int rank;
202 int retval = MPI_SUCCESS;
203 smpi_mpi_request_t srequest;
204 smpi_mpi_request_t rrequest;
205
206           rank = smpi_mpi_comm_rank(comm);
207
208           /* send */
209           /* -------------*/
210             retval = smpi_create_request(sendbuf, sendcount, sendtype, 
211                                 rank,dest,sendtag, 
212                                 comm, &srequest);
213           printf("[%d] isend request src=%d -> dst=%d (retval=%d)\n",rank,rank,dest,retval);
214           smpi_mpi_isend(srequest);
215         
216           
217           //retval = MPI_Isend( sendbuf, sendcount, sendtype, dest, sendtag, MPI_COMM_WORLD, &srequest);
218
219
220           /* recv */
221           retval = smpi_create_request(recvbuf, recvcount, recvtype, 
222                                 source, rank,recvtag, 
223                                 comm, &rrequest);
224           printf("[%d] irecv request src=%d -> dst=%d (retval=%d)\n",rank,source,rank,retval);
225           smpi_mpi_irecv(rrequest);
226
227           //retval = MPI_Irecv( recvbuf, recvcount, recvtype, source, recvtag, MPI_COMM_WORLD, &rrequest);
228
229
230           smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
231           printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
232
233           smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
234           printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
235
236           return(retval);
237 }
238
239
240 /**
241  * MPI_Wait and friends
242  **/
243 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
244 {
245   return smpi_mpi_wait(*request, status);
246 }
247
248 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
249 {
250   return smpi_mpi_waitall(count, requests, status);
251 }
252
253 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
254                      MPI_Status status[])
255 {
256   return smpi_mpi_waitany(count, requests, index, status);
257 }
258
259 /**
260  * MPI_Bcast
261  **/
262
263 /**
264  * flat bcast 
265  **/
266 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
267 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
268                 MPI_Comm comm)
269 {
270         int rank;
271         int retval = MPI_SUCCESS;
272         smpi_mpi_request_t request;
273
274         rank = smpi_mpi_comm_rank(comm);
275         if (rank == root) {
276                 retval = smpi_create_request(buf, count, datatype, root,
277                                 (root + 1) % comm->size, 0, comm, &request);
278                 request->forward = comm->size - 1;
279                 smpi_mpi_isend(request);
280         } else {
281                 retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
282                                 0, comm, &request);
283                 smpi_mpi_irecv(request);
284         }
285
286         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
287         xbt_mallocator_release(smpi_global->request_mallocator, request);
288
289         return(retval);
290
291 }
292
293 /**
294  * Bcast user entry point
295  **/
296 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
297                    MPI_Comm comm)
298 {
299   int retval = MPI_SUCCESS;
300
301   smpi_bench_end();
302
303   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
304   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
305
306   smpi_bench_begin();
307
308   return retval;
309 }
310
311
312
313 //#ifdef DEBUG_REDUCE
314 /**
315  * debugging helper function
316  **/
317 static void print_buffer_int(void *buf, int len, char *msg, int rank)
318 {
319   int tmp, *v;
320   printf("**[%d] %s: ", rank, msg);
321   for (tmp = 0; tmp < len; tmp++) {
322     v = buf;
323     printf("[%d]", v[tmp]);
324   }
325   printf("\n");
326   free(msg);
327 }
328 static void print_buffer_double(void *buf, int len, char *msg, int rank)
329 {
330   int tmp;
331   double *v;
332   printf("**[%d] %s: ", rank, msg);
333   for (tmp = 0; tmp < len; tmp++) {
334     v = buf;
335     printf("[%lf]", v[tmp]);
336   }
337   printf("\n");
338   free(msg);
339 }
340
341
342 //#endif
343 /**
344  * MPI_Reduce
345  **/
346 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
347                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
348 {
349         int retval = MPI_SUCCESS;
350         int rank;
351         int size;
352         int i;
353         int tag = 0;
354         smpi_mpi_request_t *requests;
355         smpi_mpi_request_t request;
356
357         smpi_bench_end();
358
359         rank = smpi_mpi_comm_rank(comm);
360         size = comm->size;
361
362         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
363
364 #ifdef DEBUG_REDUCE
365                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
366 #endif
367                 retval =
368                         smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
369                                         &request);
370                 smpi_mpi_isend(request);
371                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
372                 xbt_mallocator_release(smpi_global->request_mallocator, request);
373
374         } else {
375                 // i am the ROOT: wait for all buffers by creating one request by sender
376                 int src;
377                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
378
379                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
380                 for (i = 0; i < size-1; i++) {
381                         // we need 1 buffer per request to store intermediate receptions
382                         tmpbufs[i] = xbt_malloc(count * datatype->size);
383                 }  
384                 // root: initiliaze recv buf with my own snd buf
385                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
386
387                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
388                 // since we should op values as soon as one receiving request matches.
389                 for (i = 0; i < size-1; i++) {
390                         // reminder: for smpi_create_request() the src is always the process sending.
391                         src = i < root ? i : i + 1;
392                         retval = smpi_create_request(tmpbufs[i], count, datatype,
393                                         src, root, tag, comm, &(requests[i]));
394                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
395                                 if (MPI_SUCCESS == retval) {
396                                         smpi_mpi_irecv(requests[i]);
397                                 }
398                         }
399                 }
400                 // now, wait for completion of all irecv's.
401                 for (i = 0; i < size-1; i++) {
402                         int index = MPI_UNDEFINED;
403                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
404 #ifdef DEBUG_REDUCE
405                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
406                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
407                                         rank);
408 #endif
409
410                         // arg 2 is modified
411                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
412 #ifdef DEBUG_REDUCE
413                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
414 #endif
415                         xbt_free(tmpbufs[index]);
416                         /* FIXME: with the following line, it  generates an
417                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
418                          */
419                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
420                 }
421                 xbt_free(requests);
422                 xbt_free(tmpbufs);
423         }
424         smpi_bench_begin();
425         return retval;
426 }
427
428 /**
429  * MPI_Allreduce
430  *
431  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
432  **/
433 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
434                          MPI_Op op, MPI_Comm comm );
435 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
436                          MPI_Op op, MPI_Comm comm )
437 {
438   int retval = MPI_SUCCESS;
439   int root=1;  // arbitrary choice
440
441   //smpi_bench_end(); //FIXME: restaure after calling smpi_mpi_reduce instead
442
443   DEBUG0("Reduce");
444   retval = SMPI_MPI_Reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
445   if (MPI_SUCCESS != retval)
446             return(retval);
447
448   DEBUG0("Reduce done, time to bcast");
449   retval = SMPI_MPI_Bcast( sendbuf, count, datatype, root, comm);
450 //  smpi_bench_begin();
451   return( retval );
452 }
453
454
455 /**
456  * MPI_Scatter user entry point
457  **/
458 //int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
459 //                       void *recvbuf, int recvcount, MPI_Datatype recvtype,int root,
460 //                      MPI_Comm comm);
461 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
462                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
463                            int root, MPI_Comm comm)
464 {
465   int retval = MPI_SUCCESS;
466   int i;
467   int cnt=0;  
468   int rank;
469   int tag=0;
470   char *cptr;  // to manipulate the void * buffers
471   smpi_mpi_request_t *requests;
472   smpi_mpi_request_t request;
473   smpi_mpi_status_t status;
474
475
476   smpi_bench_end();
477
478   rank = smpi_mpi_comm_rank(comm);
479
480   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
481   if (rank == root) {
482           // i am the root: distribute my sendbuf
483           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
484           cptr = sendbuf;
485           for (i=0; i < comm->size; i++) {
486                   if ( i!=root ) { // send to processes ...
487
488                           retval = smpi_create_request((void *)cptr, sendcount, 
489                                           datatype, root, i, tag, comm, &(requests[cnt]));
490                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
491                                   if (MPI_SUCCESS == retval) {
492                                           smpi_mpi_isend(requests[cnt]);
493                                   }
494                                   }
495                                   cnt++;
496                         } 
497                         else { // ... except if it's me.
498                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
499                         }
500                   cptr += sendcount*datatype->size;
501             }
502             for(i=0; i<cnt; i++) { // wait for send to complete
503                             /* FIXME: waitall() should be slightly better */
504                             smpi_mpi_wait(requests[i], &status);
505                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
506
507             }
508   } 
509   else {  // i am a non-root process: wait data from the root
510             retval = smpi_create_request(recvbuf,recvcount, 
511                                   recvtype, root, rank, tag, comm, &request);
512             if (NULL != request && MPI_SUCCESS == retval) {
513                         if (MPI_SUCCESS == retval) {
514                                   smpi_mpi_irecv(request);
515                         }
516             }
517             smpi_mpi_wait(request, &status);
518             xbt_mallocator_release(smpi_global->request_mallocator, request);
519   }
520   xbt_free(requests);
521
522   smpi_bench_begin();
523
524   return retval;
525 }
526
527
528 /**
529  * MPI_Alltoall user entry point
530  * 
531  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
532  * ompi/mca/coll/tuned/coll_tuned_module.c
533  **/
534 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
535                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
536                            MPI_Comm comm)
537 {
538   int retval = MPI_SUCCESS;
539   int block_dsize;
540   int rank;
541
542   smpi_bench_end();
543
544   rank = smpi_mpi_comm_rank(comm);
545   block_dsize = datatype->size * sendcount;
546
547   if ((block_dsize < 200) && (comm->size > 12)) {
548             retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
549                                   recvbuf, recvcount, recvtype, comm);
550
551   } else if (block_dsize < 3000) {
552 /* use this one !!          retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, datatype,
553                                   recvbuf, recvcount, recvtype, comm);
554                                   */
555   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
556                                   recvbuf, recvcount, recvtype, comm);
557   } else {
558
559   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
560                                   recvbuf, recvcount, recvtype, comm);
561   }
562
563   smpi_bench_begin();
564
565   return retval;
566 }
567
568
569
570
571 // used by comm_split to sort ranks based on key values
572 int smpi_compare_rankkeys(const void *a, const void *b);
573 int smpi_compare_rankkeys(const void *a, const void *b)
574 {
575   int *x = (int *) a;
576   int *y = (int *) b;
577
578   if (x[1] < y[1])
579     return -1;
580
581   if (x[1] == y[1]) {
582     if (x[0] < y[0])
583       return -1;
584     if (x[0] == y[0])
585       return 0;
586     return 1;
587   }
588
589   return 1;
590 }
591
592 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
593                         MPI_Comm * comm_out)
594 {
595   int retval = MPI_SUCCESS;
596
597   int index, rank;
598   smpi_mpi_request_t request;
599   int colorkey[2];
600   smpi_mpi_status_t status;
601
602   smpi_bench_end();
603
604   // FIXME: need to test parameters
605
606   index = smpi_process_index();
607   rank = comm->index_to_rank_map[index];
608
609   // default output
610   comm_out = NULL;
611
612   // root node does most of the real work
613   if (0 == rank) {
614     int colormap[comm->size];
615     int keymap[comm->size];
616     int rankkeymap[comm->size * 2];
617     int i, j;
618     smpi_mpi_communicator_t tempcomm = NULL;
619     int count;
620     int indextmp;
621
622     colormap[0] = color;
623     keymap[0] = key;
624
625     // FIXME: use scatter/gather or similar instead of individual comms
626     for (i = 1; i < comm->size; i++) {
627       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
628                                    rank, MPI_ANY_TAG, comm, &request);
629       smpi_mpi_irecv(request);
630       smpi_mpi_wait(request, &status);
631       colormap[status.MPI_SOURCE] = colorkey[0];
632       keymap[status.MPI_SOURCE] = colorkey[1];
633       xbt_mallocator_release(smpi_global->request_mallocator, request);
634     }
635
636     for (i = 0; i < comm->size; i++) {
637       if (MPI_UNDEFINED == colormap[i]) {
638         continue;
639       }
640       // make a list of nodes with current color and sort by keys
641       count = 0;
642       for (j = i; j < comm->size; j++) {
643         if (colormap[i] == colormap[j]) {
644           colormap[j] = MPI_UNDEFINED;
645           rankkeymap[count * 2] = j;
646           rankkeymap[count * 2 + 1] = keymap[j];
647           count++;
648         }
649       }
650       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
651
652       // new communicator
653       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
654       tempcomm->barrier_count = 0;
655       tempcomm->size = count;
656       tempcomm->barrier_mutex = SIMIX_mutex_init();
657       tempcomm->barrier_cond = SIMIX_cond_init();
658       tempcomm->rank_to_index_map = xbt_new(int, count);
659       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
660       for (j = 0; j < smpi_global->process_count; j++) {
661         tempcomm->index_to_rank_map[j] = -1;
662       }
663       for (j = 0; j < count; j++) {
664         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
665         tempcomm->rank_to_index_map[j] = indextmp;
666         tempcomm->index_to_rank_map[indextmp] = j;
667       }
668       for (j = 0; j < count; j++) {
669         if (rankkeymap[j * 2]) {
670           retval = smpi_create_request(&j, 1, MPI_INT, 0,
671                                        rankkeymap[j * 2], 0, comm, &request);
672           request->data = tempcomm;
673           smpi_mpi_isend(request);
674           smpi_mpi_wait(request, &status);
675           xbt_mallocator_release(smpi_global->request_mallocator, request);
676         } else {
677           *comm_out = tempcomm;
678         }
679       }
680     }
681   } else {
682     colorkey[0] = color;
683     colorkey[1] = key;
684     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
685                                  &request);
686     smpi_mpi_isend(request);
687     smpi_mpi_wait(request, &status);
688     xbt_mallocator_release(smpi_global->request_mallocator, request);
689     if (MPI_UNDEFINED != color) {
690       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
691                                    &request);
692       smpi_mpi_irecv(request);
693       smpi_mpi_wait(request, &status);
694       *comm_out = request->data;
695     }
696   }
697
698   smpi_bench_begin();
699
700   return retval;
701 }
702
703 double SMPI_MPI_Wtime(void)
704 {
705   return (SIMIX_get_clock());
706 }
707
708