Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
08d3eba5f4acc3a82e6110a1a655a443dbf511cd
[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 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
190                    MPI_Comm comm)
191 {
192
193   int retval = MPI_SUCCESS;
194   int rank;
195   smpi_mpi_request_t request;
196
197   smpi_bench_end();
198
199   rank = smpi_mpi_comm_rank(comm);
200
201   if (rank == root) {
202     retval = smpi_create_request(buf, count, datatype, root,
203                                  (root + 1) % comm->size, 0, comm, &request);
204     request->forward = comm->size - 1;
205     smpi_mpi_isend(request);
206   } else {
207     retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
208                                  0, comm, &request);
209     smpi_mpi_irecv(request);
210   }
211
212   smpi_mpi_wait(request, MPI_STATUS_IGNORE);
213   xbt_mallocator_release(smpi_global->request_mallocator, request);
214
215   smpi_bench_begin();
216
217   return retval;
218 }
219
220 // used by comm_split to sort ranks based on key values
221 int smpi_compare_rankkeys(const void *a, const void *b);
222 int smpi_compare_rankkeys(const void *a, const void *b)
223 {
224   int *x = (int *) a;
225   int *y = (int *) b;
226
227   if (x[1] < y[1])
228     return -1;
229
230   if (x[1] == y[1]) {
231     if (x[0] < y[0])
232       return -1;
233     if (x[0] == y[0])
234       return 0;
235     return 1;
236   }
237
238   return 1;
239 }
240
241 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
242                         MPI_Comm * comm_out)
243 {
244   int retval = MPI_SUCCESS;
245
246   int index, rank;
247   smpi_mpi_request_t request;
248   int colorkey[2];
249   smpi_mpi_status_t status;
250
251   smpi_bench_end();
252
253   // FIXME: need to test parameters
254
255   index = smpi_process_index();
256   rank = comm->index_to_rank_map[index];
257
258   // default output
259   comm_out = NULL;
260
261   // root node does most of the real work
262   if (0 == rank) {
263     int colormap[comm->size];
264     int keymap[comm->size];
265     int rankkeymap[comm->size * 2];
266     int i, j;
267     smpi_mpi_communicator_t tempcomm = NULL;
268     int count;
269     int indextmp;
270
271     colormap[0] = color;
272     keymap[0] = key;
273
274     // FIXME: use scatter/gather or similar instead of individual comms
275     for (i = 1; i < comm->size; i++) {
276       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
277                                    rank, MPI_ANY_TAG, comm, &request);
278       smpi_mpi_irecv(request);
279       smpi_mpi_wait(request, &status);
280       colormap[status.MPI_SOURCE] = colorkey[0];
281       keymap[status.MPI_SOURCE] = colorkey[1];
282       xbt_mallocator_release(smpi_global->request_mallocator, request);
283     }
284
285     for (i = 0; i < comm->size; i++) {
286       if (MPI_UNDEFINED == colormap[i]) {
287         continue;
288       }
289       // make a list of nodes with current color and sort by keys
290       count = 0;
291       for (j = i; j < comm->size; j++) {
292         if (colormap[i] == colormap[j]) {
293           colormap[j] = MPI_UNDEFINED;
294           rankkeymap[count * 2] = j;
295           rankkeymap[count * 2 + 1] = keymap[j];
296           count++;
297         }
298       }
299       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
300
301       // new communicator
302       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
303       tempcomm->barrier_count = 0;
304       tempcomm->size = count;
305       tempcomm->barrier_mutex = SIMIX_mutex_init();
306       tempcomm->barrier_cond = SIMIX_cond_init();
307       tempcomm->rank_to_index_map = xbt_new(int, count);
308       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
309       for (j = 0; j < smpi_global->process_count; j++) {
310         tempcomm->index_to_rank_map[j] = -1;
311       }
312       for (j = 0; j < count; j++) {
313         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
314         tempcomm->rank_to_index_map[j] = indextmp;
315         tempcomm->index_to_rank_map[indextmp] = j;
316       }
317       for (j = 0; j < count; j++) {
318         if (rankkeymap[j * 2]) {
319           retval = smpi_create_request(&j, 1, MPI_INT, 0,
320                                        rankkeymap[j * 2], 0, comm, &request);
321           request->data = tempcomm;
322           smpi_mpi_isend(request);
323           smpi_mpi_wait(request, &status);
324           xbt_mallocator_release(smpi_global->request_mallocator, request);
325         } else {
326           *comm_out = tempcomm;
327         }
328       }
329     }
330   } else {
331     colorkey[0] = color;
332     colorkey[1] = key;
333     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
334                                  &request);
335     smpi_mpi_isend(request);
336     smpi_mpi_wait(request, &status);
337     xbt_mallocator_release(smpi_global->request_mallocator, request);
338     if (MPI_UNDEFINED != color) {
339       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
340                                    &request);
341       smpi_mpi_irecv(request);
342       smpi_mpi_wait(request, &status);
343       *comm_out = request->data;
344     }
345   }
346
347   smpi_bench_begin();
348
349   return retval;
350 }
351
352 double SMPI_MPI_Wtime( void ) 
353
354           return ( SIMIX_get_clock() );
355 }