Logo AND Algorithmique Numérique Distribuée

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