Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Make Waitall and Waitany visible from user side (and fix their prototype to stick...
[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 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
185 {
186   return smpi_mpi_wait(*request, status);
187 }
188
189 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[]) {
190         return smpi_mpi_waitall(count, requests,status);
191 }
192 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status status[]) {
193         return smpi_mpi_waitany(count, requests, index,status);
194 }
195 /**
196  * MPI_Bcast
197  **/
198 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
199                    MPI_Comm comm)
200 {
201
202   int retval = MPI_SUCCESS;
203   int rank;
204   smpi_mpi_request_t request;
205
206   smpi_bench_end();
207
208   rank = smpi_mpi_comm_rank(comm);
209
210   if (rank == root) {
211     retval = smpi_create_request(buf, count, datatype, root,
212                                  (root + 1) % comm->size, 0, comm, &request);
213     request->forward = comm->size - 1;
214     smpi_mpi_isend(request);
215   } else {
216     retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
217                                  0, comm, &request);
218     smpi_mpi_irecv(request);
219   }
220
221   smpi_mpi_wait(request, MPI_STATUS_IGNORE);
222   xbt_mallocator_release(smpi_global->request_mallocator, request);
223
224   smpi_bench_begin();
225
226   return retval;
227 }
228
229 /**
230  * MPI_Reduce
231  **/
232
233 int SMPI_MPI_Reduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
234                          MPI_Op op, int root, MPI_Comm comm )
235 {
236   int retval = MPI_SUCCESS;
237   int rank;
238   int size;
239   int i;
240   smpi_mpi_request_t *tabrequest;
241
242   smpi_bench_end();
243
244   rank = smpi_mpi_comm_rank(comm);
245   size = comm->size;
246
247   tabrequest = xbt_malloc((size)*sizeof(smpi_mpi_request_t));
248
249   if (rank != root) { // if i am not root, simply send my buffer to root
250             retval = smpi_create_request(sendbuf, count, datatype,
251                                   rank, root, 0, comm, &(tabrequest[rank]));
252             smpi_mpi_isend(tabrequest[rank]);
253             smpi_mpi_wait(tabrequest[rank], MPI_STATUS_IGNORE);
254             //printf("DEBUG: rank %d sent my sendbuf to root (rank %d)\n",rank,root);
255   } else {
256             // i am the root: wait for all buffers by creating requests
257             // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
258             // since we should op values as soon as one receiving request matches.
259             for (i=0; i<comm->size; i++) {
260                         if ( rank != i ) { // except for me
261                                   // reminder: for smpi_create_request() the src is always the process sending.
262                                   retval = smpi_create_request(recvbuf, count, datatype, MPI_ANY_SOURCE, root,
263                                                         0, comm, &(tabrequest[i]));
264                                   if (NULL != tabrequest[i] && MPI_SUCCESS == retval) {
265                                             if (MPI_SUCCESS == retval) {
266                                                         smpi_mpi_irecv(tabrequest[i]);
267                                             }
268                                   }
269                         }
270             }
271             // now, wait for completion of all irecv's.
272             // FIXME: we should implement smpi_wait_all for a more asynchronous behavior
273             for (i=0; i<comm->size; i++) {
274                         if ( rank != i ) { // except for me
275                                   smpi_mpi_wait(tabrequest[i], MPI_STATUS_IGNORE);
276
277                                   // FIXME: the core part is here. To be written ...
278
279                                   fprintf(stderr,"[smpi] %s:%d : MPI_Reduce *Not yet implemented*.\n",__FILE__,__LINE__);
280                         }
281             }
282
283   }
284   for (i=0; i<comm->size; i++)
285             xbt_mallocator_release(smpi_global->request_mallocator, tabrequest[i]);
286
287   smpi_bench_begin();
288
289   return retval;
290 }
291
292 // used by comm_split to sort ranks based on key values
293 int smpi_compare_rankkeys(const void *a, const void *b);
294 int smpi_compare_rankkeys(const void *a, const void *b)
295 {
296   int *x = (int *) a;
297   int *y = (int *) b;
298
299   if (x[1] < y[1])
300     return -1;
301
302   if (x[1] == y[1]) {
303     if (x[0] < y[0])
304       return -1;
305     if (x[0] == y[0])
306       return 0;
307     return 1;
308   }
309
310   return 1;
311 }
312
313 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
314                         MPI_Comm * comm_out)
315 {
316   int retval = MPI_SUCCESS;
317
318   int index, rank;
319   smpi_mpi_request_t request;
320   int colorkey[2];
321   smpi_mpi_status_t status;
322
323   smpi_bench_end();
324
325   // FIXME: need to test parameters
326
327   index = smpi_process_index();
328   rank = comm->index_to_rank_map[index];
329
330   // default output
331   comm_out = NULL;
332
333   // root node does most of the real work
334   if (0 == rank) {
335     int colormap[comm->size];
336     int keymap[comm->size];
337     int rankkeymap[comm->size * 2];
338     int i, j;
339     smpi_mpi_communicator_t tempcomm = NULL;
340     int count;
341     int indextmp;
342
343     colormap[0] = color;
344     keymap[0] = key;
345
346     // FIXME: use scatter/gather or similar instead of individual comms
347     for (i = 1; i < comm->size; i++) {
348       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
349                                    rank, MPI_ANY_TAG, comm, &request);
350       smpi_mpi_irecv(request);
351       smpi_mpi_wait(request, &status);
352       colormap[status.MPI_SOURCE] = colorkey[0];
353       keymap[status.MPI_SOURCE] = colorkey[1];
354       xbt_mallocator_release(smpi_global->request_mallocator, request);
355     }
356
357     for (i = 0; i < comm->size; i++) {
358       if (MPI_UNDEFINED == colormap[i]) {
359         continue;
360       }
361       // make a list of nodes with current color and sort by keys
362       count = 0;
363       for (j = i; j < comm->size; j++) {
364         if (colormap[i] == colormap[j]) {
365           colormap[j] = MPI_UNDEFINED;
366           rankkeymap[count * 2] = j;
367           rankkeymap[count * 2 + 1] = keymap[j];
368           count++;
369         }
370       }
371       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
372
373       // new communicator
374       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
375       tempcomm->barrier_count = 0;
376       tempcomm->size = count;
377       tempcomm->barrier_mutex = SIMIX_mutex_init();
378       tempcomm->barrier_cond = SIMIX_cond_init();
379       tempcomm->rank_to_index_map = xbt_new(int, count);
380       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
381       for (j = 0; j < smpi_global->process_count; j++) {
382         tempcomm->index_to_rank_map[j] = -1;
383       }
384       for (j = 0; j < count; j++) {
385         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
386         tempcomm->rank_to_index_map[j] = indextmp;
387         tempcomm->index_to_rank_map[indextmp] = j;
388       }
389       for (j = 0; j < count; j++) {
390         if (rankkeymap[j * 2]) {
391           retval = smpi_create_request(&j, 1, MPI_INT, 0,
392                                        rankkeymap[j * 2], 0, comm, &request);
393           request->data = tempcomm;
394           smpi_mpi_isend(request);
395           smpi_mpi_wait(request, &status);
396           xbt_mallocator_release(smpi_global->request_mallocator, request);
397         } else {
398           *comm_out = tempcomm;
399         }
400       }
401     }
402   } else {
403     colorkey[0] = color;
404     colorkey[1] = key;
405     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
406                                  &request);
407     smpi_mpi_isend(request);
408     smpi_mpi_wait(request, &status);
409     xbt_mallocator_release(smpi_global->request_mallocator, request);
410     if (MPI_UNDEFINED != color) {
411       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
412                                    &request);
413       smpi_mpi_irecv(request);
414       smpi_mpi_wait(request, &status);
415       *comm_out = request->data;
416     }
417   }
418
419   smpi_bench_begin();
420
421   return retval;
422 }
423
424 double SMPI_MPI_Wtime( void )
425 {
426           return ( SIMIX_get_clock() );
427 }