Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
* minimum of datatype handling for alltoall tuned
[simgrid.git] / src / smpi / smpi_mpi.c
1
2
3 #include "private.h"
4 #include "smpi_coll_private.h"
5 #include "smpi_mpi_dt_private.h"
6
7 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
8                                 "Logging specific to SMPI (mpi)");
9
10 int SMPI_MPI_Init(int *argc, char ***argv)
11 {
12   smpi_process_init(argc, argv);
13   smpi_bench_begin();
14   return MPI_SUCCESS;
15 }
16
17 int SMPI_MPI_Finalize()
18 {
19   smpi_bench_end();
20   smpi_process_finalize();
21   return MPI_SUCCESS;
22 }
23
24 // right now this just exits the current node, should send abort signal to all
25 // hosts in the communicator (TODO)
26 int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
27 {
28   smpi_exit(errorcode);
29   return 0;
30 }
31
32 int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
33 {
34   int retval = MPI_SUCCESS;
35
36   smpi_bench_end();
37
38   if (NULL == comm) {
39     retval = MPI_ERR_COMM;
40   } else if (NULL == size) {
41     retval = MPI_ERR_ARG;
42   } else {
43     *size = comm->size;
44   }
45
46   smpi_bench_begin();
47
48   return retval;
49 }
50
51 int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
52 {
53   int retval = MPI_SUCCESS;
54
55   smpi_bench_end();
56
57   if (NULL == comm) {
58     retval = MPI_ERR_COMM;
59   } else if (NULL == rank) {
60     retval = MPI_ERR_ARG;
61   } else {
62     *rank = smpi_mpi_comm_rank(comm);
63   }
64
65   smpi_bench_begin();
66
67   return retval;
68 }
69
70
71
72 /**
73  * Barrier
74  **/
75 int SMPI_MPI_Barrier(MPI_Comm comm)
76 {
77   int retval = MPI_SUCCESS;
78   int arity=4;
79
80   smpi_bench_end();
81
82   if (NULL == comm) {
83     retval = MPI_ERR_COMM;
84   } else {
85
86     /*
87      * original implemantation:
88      * retval = smpi_mpi_barrier(comm);
89      * this one is unrealistic: it just cond_waits, means no time.
90      */
91      retval = nary_tree_barrier( comm, arity );
92   }
93
94   smpi_bench_begin();
95
96   return retval;
97 }
98
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 /**
161  * MPI_Send user level
162  **/
163 int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
164                   int tag, MPI_Comm comm)
165 {
166   int retval = MPI_SUCCESS;
167   smpi_mpi_request_t request;
168
169   smpi_bench_end();
170
171   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
172                                &request);
173   if (NULL != request && MPI_SUCCESS == retval) {
174     retval = smpi_mpi_isend(request);
175     if (MPI_SUCCESS == retval) {
176       smpi_mpi_wait(request, MPI_STATUS_IGNORE);
177     }
178     xbt_mallocator_release(smpi_global->request_mallocator, request);
179   }
180
181   smpi_bench_begin();
182
183   return retval;
184 }
185
186
187 /**
188  * MPI_Sendrecv internal level 
189  **/
190 int smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
191                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
192                     MPI_Comm comm, MPI_Status *status)
193 {
194 int rank;
195 int retval = MPI_SUCCESS;
196 smpi_mpi_request_t srequest;
197 smpi_mpi_request_t rrequest;
198
199           rank = smpi_mpi_comm_rank(comm);
200
201           /* send */
202           retval = smpi_create_request(sendbuf, sendcount, sendtype, 
203                                 rank,dest,sendtag, 
204                                 comm, &srequest);
205           smpi_mpi_isend(srequest);
206         
207           /* recv */
208           retval = smpi_create_request(recvbuf, recvcount, recvtype, 
209                                 source, rank,recvtag, 
210                                 comm, &rrequest);
211           smpi_mpi_irecv(rrequest);
212
213           smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
214           //printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
215           smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
216           //printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
217
218           return(retval);
219 }
220 /**
221  * MPI_Sendrecv user entry point
222  **/
223 int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
224                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
225                     MPI_Comm comm, MPI_Status *status)
226 {
227 int retval = MPI_SUCCESS;
228
229   smpi_bench_end();
230   smpi_mpi_sendrecv( sendbuf, sendcount, sendtype, dest, sendtag, 
231                          recvbuf, recvcount, recvtype, source, recvtag,
232                          comm, status);
233   smpi_bench_begin();
234
235   return retval;
236
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  * Bcast internal level 
294  **/
295 int smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
296                    MPI_Comm comm)
297 {
298   int retval = MPI_SUCCESS;
299   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
300   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
301   return retval;
302 }
303
304 /**
305  * Bcast user entry point
306  **/
307 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
308                    MPI_Comm comm)
309 {
310   int retval = MPI_SUCCESS;
311
312   smpi_bench_end();
313   smpi_mpi_bcast(buf,count,datatype,root,comm);
314   smpi_bench_begin();
315
316   return retval;
317 }
318
319
320
321 #ifdef DEBUG_REDUCE
322 /**
323  * debugging helper function
324  **/
325 static void print_buffer_int(void *buf, int len, char *msg, int rank)
326 {
327   int tmp, *v;
328   printf("**[%d] %s: ", rank, msg);
329   for (tmp = 0; tmp < len; tmp++) {
330     v = buf;
331     printf("[%d]", v[tmp]);
332   }
333   printf("\n");
334   free(msg);
335 }
336 static void print_buffer_double(void *buf, int len, char *msg, int rank)
337 {
338   int tmp;
339   double *v;
340   printf("**[%d] %s: ", rank, msg);
341   for (tmp = 0; tmp < len; tmp++) {
342     v = buf;
343     printf("[%lf]", v[tmp]);
344   }
345   printf("\n");
346   free(msg);
347 }
348
349
350 #endif
351 /**
352  * MPI_Reduce internal level 
353  **/
354 int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
355                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
356 {
357         int retval = MPI_SUCCESS;
358         int rank;
359         int size;
360         int i;
361         int tag = 0;
362         smpi_mpi_request_t *requests;
363         smpi_mpi_request_t request;
364
365         smpi_bench_end();
366
367         rank = smpi_mpi_comm_rank(comm);
368         size = comm->size;
369
370         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
371
372 #ifdef DEBUG_REDUCE
373                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
374 #endif
375                 retval = smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
376                                         &request);
377                 smpi_mpi_isend(request);
378                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
379                 xbt_mallocator_release(smpi_global->request_mallocator, request);
380
381         } else {
382                 // i am the ROOT: wait for all buffers by creating one request by sender
383                 int src;
384                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
385
386                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
387                 for (i = 0; i < size-1; i++) {
388                         // we need 1 buffer per request to store intermediate receptions
389                         tmpbufs[i] = xbt_malloc(count * datatype->size);
390                 }  
391                 // root: initiliaze recv buf with my own snd buf
392                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
393
394                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
395                 // since we should op values as soon as one receiving request matches.
396                 for (i = 0; i < size-1; i++) {
397                         // reminder: for smpi_create_request() the src is always the process sending.
398                         src = i < root ? i : i + 1;
399                         retval = smpi_create_request(tmpbufs[i], count, datatype,
400                                         src, root, tag, comm, &(requests[i]));
401                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
402                                 if (MPI_SUCCESS == retval) {
403                                         smpi_mpi_irecv(requests[i]);
404                                 }
405                         }
406                 }
407                 // now, wait for completion of all irecv's.
408                 for (i = 0; i < size-1; i++) {
409                         int index = MPI_UNDEFINED;
410                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
411 #ifdef DEBUG_REDUCE
412                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
413                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
414                                         rank);
415 #endif
416
417                         // arg 2 is modified
418                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
419 #ifdef DEBUG_REDUCE
420                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
421 #endif
422                         xbt_free(tmpbufs[index]);
423                         /* FIXME: with the following line, it  generates an
424                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
425                          */
426                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
427                 }
428                 xbt_free(requests);
429                 xbt_free(tmpbufs);
430         }
431         return retval;
432 }
433
434 /**
435  * MPI_Reduce user entry point
436  **/
437 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
438                     MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
439 {
440 int retval = MPI_SUCCESS;
441
442           smpi_bench_end();
443
444           retval = smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
445
446           smpi_bench_begin();
447           return retval;
448 }
449
450
451
452 /**
453  * MPI_Allreduce
454  *
455  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
456  **/
457 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
458                     MPI_Op op, MPI_Comm comm )
459 {
460 int retval = MPI_SUCCESS;
461 int root=0;  // arbitrary choice
462
463           smpi_bench_end();
464
465           retval = smpi_mpi_reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
466           if (MPI_SUCCESS != retval)
467                     return(retval);
468
469           retval = smpi_mpi_bcast( sendbuf, count, datatype, root, comm);
470
471           smpi_bench_end();
472           return( retval );
473 }
474
475
476 /**
477  * MPI_Scatter user entry point
478  **/
479 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
480                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
481                            int root, MPI_Comm comm)
482 {
483   int retval = MPI_SUCCESS;
484   int i;
485   int cnt=0;  
486   int rank;
487   int tag=0;
488   char *cptr;  // to manipulate the void * buffers
489   smpi_mpi_request_t *requests;
490   smpi_mpi_request_t request;
491   smpi_mpi_status_t status;
492
493
494   smpi_bench_end();
495
496   rank = smpi_mpi_comm_rank(comm);
497
498   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
499   if (rank == root) {
500           // i am the root: distribute my sendbuf
501           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
502           cptr = sendbuf;
503           for (i=0; i < comm->size; i++) {
504                   if ( i!=root ) { // send to processes ...
505
506                           retval = smpi_create_request((void *)cptr, sendcount, 
507                                           datatype, root, i, tag, comm, &(requests[cnt]));
508                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
509                                   if (MPI_SUCCESS == retval) {
510                                           smpi_mpi_isend(requests[cnt]);
511                                   }
512                                   }
513                                   cnt++;
514                         } 
515                         else { // ... except if it's me.
516                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
517                         }
518                   cptr += sendcount*datatype->size;
519             }
520             for(i=0; i<cnt; i++) { // wait for send to complete
521                             /* FIXME: waitall() should be slightly better */
522                             smpi_mpi_wait(requests[i], &status);
523                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
524
525             }
526   } 
527   else {  // i am a non-root process: wait data from the root
528             retval = smpi_create_request(recvbuf,recvcount, 
529                                   recvtype, root, rank, tag, comm, &request);
530             if (NULL != request && MPI_SUCCESS == retval) {
531                         if (MPI_SUCCESS == retval) {
532                                   smpi_mpi_irecv(request);
533                         }
534             }
535             smpi_mpi_wait(request, &status);
536             xbt_mallocator_release(smpi_global->request_mallocator, request);
537   }
538   xbt_free(requests);
539
540   smpi_bench_begin();
541
542   return retval;
543 }
544
545
546 /**
547  * MPI_Alltoall user entry point
548  * 
549  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
550  * ompi/mca/coll/tuned/coll_tuned_module.c
551  **/
552 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
553                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
554                            MPI_Comm comm)
555 {
556   int retval = MPI_SUCCESS;
557   int block_dsize;
558   int rank;
559
560   smpi_bench_end();
561
562   rank = smpi_mpi_comm_rank(comm);
563   block_dsize = datatype->size * sendcount;
564   INFO2("[%d] optimized alltoall() called. Block size sent to each rank=%d.\n",rank,block_dsize);
565
566   if ((block_dsize < 200) && (comm->size > 12)) {
567             retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
568                                   recvbuf, recvcount, recvtype, comm);
569
570   } else if (block_dsize < 3000) {
571             retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, datatype,
572                                   recvbuf, recvcount, recvtype, comm);
573   } else {
574
575   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
576                                   recvbuf, recvcount, recvtype, comm);
577   }
578
579   smpi_bench_begin();
580
581   return retval;
582 }
583
584
585
586
587 // used by comm_split to sort ranks based on key values
588 int smpi_compare_rankkeys(const void *a, const void *b);
589 int smpi_compare_rankkeys(const void *a, const void *b)
590 {
591   int *x = (int *) a;
592   int *y = (int *) b;
593
594   if (x[1] < y[1])
595     return -1;
596
597   if (x[1] == y[1]) {
598     if (x[0] < y[0])
599       return -1;
600     if (x[0] == y[0])
601       return 0;
602     return 1;
603   }
604
605   return 1;
606 }
607
608 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
609                         MPI_Comm * comm_out)
610 {
611   int retval = MPI_SUCCESS;
612
613   int index, rank;
614   smpi_mpi_request_t request;
615   int colorkey[2];
616   smpi_mpi_status_t status;
617
618   smpi_bench_end();
619
620   // FIXME: need to test parameters
621
622   index = smpi_process_index();
623   rank = comm->index_to_rank_map[index];
624
625   // default output
626   comm_out = NULL;
627
628   // root node does most of the real work
629   if (0 == rank) {
630     int colormap[comm->size];
631     int keymap[comm->size];
632     int rankkeymap[comm->size * 2];
633     int i, j;
634     smpi_mpi_communicator_t tempcomm = NULL;
635     int count;
636     int indextmp;
637
638     colormap[0] = color;
639     keymap[0] = key;
640
641     // FIXME: use scatter/gather or similar instead of individual comms
642     for (i = 1; i < comm->size; i++) {
643       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
644                                    rank, MPI_ANY_TAG, comm, &request);
645       smpi_mpi_irecv(request);
646       smpi_mpi_wait(request, &status);
647       colormap[status.MPI_SOURCE] = colorkey[0];
648       keymap[status.MPI_SOURCE] = colorkey[1];
649       xbt_mallocator_release(smpi_global->request_mallocator, request);
650     }
651
652     for (i = 0; i < comm->size; i++) {
653       if (MPI_UNDEFINED == colormap[i]) {
654         continue;
655       }
656       // make a list of nodes with current color and sort by keys
657       count = 0;
658       for (j = i; j < comm->size; j++) {
659         if (colormap[i] == colormap[j]) {
660           colormap[j] = MPI_UNDEFINED;
661           rankkeymap[count * 2] = j;
662           rankkeymap[count * 2 + 1] = keymap[j];
663           count++;
664         }
665       }
666       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
667
668       // new communicator
669       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
670       tempcomm->barrier_count = 0;
671       tempcomm->size = count;
672       tempcomm->barrier_mutex = SIMIX_mutex_init();
673       tempcomm->barrier_cond = SIMIX_cond_init();
674       tempcomm->rank_to_index_map = xbt_new(int, count);
675       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
676       for (j = 0; j < smpi_global->process_count; j++) {
677         tempcomm->index_to_rank_map[j] = -1;
678       }
679       for (j = 0; j < count; j++) {
680         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
681         tempcomm->rank_to_index_map[j] = indextmp;
682         tempcomm->index_to_rank_map[indextmp] = j;
683       }
684       for (j = 0; j < count; j++) {
685         if (rankkeymap[j * 2]) {
686           retval = smpi_create_request(&j, 1, MPI_INT, 0,
687                                        rankkeymap[j * 2], 0, comm, &request);
688           request->data = tempcomm;
689           smpi_mpi_isend(request);
690           smpi_mpi_wait(request, &status);
691           xbt_mallocator_release(smpi_global->request_mallocator, request);
692         } else {
693           *comm_out = tempcomm;
694         }
695       }
696     }
697   } else {
698     colorkey[0] = color;
699     colorkey[1] = key;
700     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
701                                  &request);
702     smpi_mpi_isend(request);
703     smpi_mpi_wait(request, &status);
704     xbt_mallocator_release(smpi_global->request_mallocator, request);
705     if (MPI_UNDEFINED != color) {
706       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
707                                    &request);
708       smpi_mpi_irecv(request);
709       smpi_mpi_wait(request, &status);
710       *comm_out = request->data;
711     }
712   }
713
714   smpi_bench_begin();
715
716   return retval;
717 }
718
719 double SMPI_MPI_Wtime(void)
720 {
721   return (SIMIX_get_clock());
722 }
723
724