Logo AND Algorithmique Numérique Distribuée

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