Logo AND Algorithmique Numérique Distribuée

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