Logo AND Algorithmique Numérique Distribuée

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