Logo AND Algorithmique Numérique Distribuée

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