Logo AND Algorithmique Numérique Distribuée

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