Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
A lot has been moved over to simix, but getting rid of tasks requires
[simgrid.git] / src / smpi / src / smpi_mpi.c
1 #include <stdio.h>
2 #include <sys/time.h>
3 #include "msg/msg.h"
4 #include "xbt/sysdep.h"
5 #include "xbt/xbt_portability.h"
6 #include "smpi.h"
7
8 int MPI_Init(int *argc, char ***argv) {
9   smpi_mpi_init();
10   smpi_bench_begin();
11   return MPI_SUCCESS;
12 }
13
14 int MPI_Finalize() {
15   smpi_bench_end();
16   smpi_mpi_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;
22 int MPI_Abort(MPI_Comm comm, int errorcode) {
23   smpi_exit(errorcode);
24 }
25
26 int MPI_Comm_size(MPI_Comm comm, int *size) {
27   int retval = MPI_SUCCESS;
28   smpi_bench_end();
29   if (NULL == comm) {
30     retval = MPI_ERR_COMM;
31   } else if (NULL == size) {
32     retval = MPI_ERR_ARG;
33   } else {
34     *size = comm->size;
35   }
36   smpi_bench_begin();
37   return retval;
38 }
39
40 int MPI_Comm_rank(MPI_Comm comm, int *rank) {
41   int retval = MPI_SUCCESS;
42   smpi_bench_end();
43   if (NULL == comm) {
44     retval = MPI_ERR_COMM;
45   } else if (NULL == rank) {
46     retval = MPI_ERR_ARG;
47   } else {
48     *rank = smpi_comm_rank(comm, SIMIX_host_self());
49   }
50   smpi_bench_begin();
51   return retval;
52 }
53
54 /*
55 int MPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm *newcomm) {
56   int retval = MPI_SUCCESS;
57   m_host_t host = SIMIX_host_self();
58   int rank = smpi_comm_rank(comm, host);
59   smpi_mpi_comm_split_table_node_t *split_table; 
60   split_table = xbt_malloc(sizeof(smpi_mpi_comm_split_table_node_t) * comm->size);
61   split_table[rank].color = color;
62   split_table[rank].key   = key;
63   split_table[rank].host  = host;
64   smpi_mpi_alltoall(
65   return retval;
66 }
67 */
68
69 int MPI_Type_size(MPI_Datatype datatype, size_t *size) {
70   int retval = MPI_SUCCESS;
71   smpi_bench_end();
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   smpi_bench_begin();
80   return retval;
81 }
82
83 int MPI_Wait(MPI_Request *request, MPI_Status *status) {
84   int retval = MPI_SUCCESS;
85   smpi_bench_end();
86   if (NULL == request) {
87     retval = MPI_ERR_REQUEST;
88   } else if (NULL == status) {
89     retval = MPI_ERR_ARG;
90   } else {
91     smpi_wait(*request, status);
92   }
93   smpi_bench_begin();
94   return retval;
95 }
96
97 int MPI_Waitall(int count, MPI_Request *requests, MPI_Status *statuses) {
98   int retval = MPI_SUCCESS;
99   smpi_bench_end();
100   if (NULL == requests) {
101     retval = MPI_ERR_REQUEST;
102   } else if (NULL == statuses) {
103     retval = MPI_ERR_ARG;
104   } else {
105     smpi_wait_all(count, requests, statuses);
106   }
107   smpi_bench_begin();
108   return MPI_ERR_INTERN;
109 }
110
111 int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request *request) {
112   int retval = MPI_SUCCESS;
113   int dst;
114   smpi_mpi_request_t *recvreq;
115   smpi_bench_end();
116   dst = smpi_comm_rank(comm, SIMIX_host_self());
117   retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &recvreq);
118   if (NULL != recvreq) {
119     smpi_irecv(recvreq);
120     if (NULL != request) {
121       *request = recvreq;
122     }
123   }
124   smpi_bench_begin();
125   return retval;
126 }
127
128 int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status *status) {
129   int retval = MPI_SUCCESS;
130   int dst;
131   smpi_mpi_request_t *recvreq;
132   smpi_bench_end();
133   dst = smpi_comm_rank(comm, SIMIX_host_self());
134   retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &recvreq);
135   if (NULL != recvreq) {
136     smpi_irecv(recvreq);
137     smpi_wait(recvreq, status);
138     xbt_free(recvreq);
139   }
140   smpi_bench_begin();
141   return retval;
142 }
143
144 int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request *request) {
145   int retval = MPI_SUCCESS;
146   int src;
147   smpi_mpi_request_t *sendreq;
148   smpi_bench_end();
149   src = smpi_comm_rank(comm, SIMIX_host_self());
150   retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &sendreq);
151   if (NULL != sendreq) {
152     smpi_isend(sendreq);
153     if (NULL != request) {
154       *request = sendreq;
155     }
156   }
157   smpi_bench_begin();
158   return retval;
159 }
160
161 int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
162   int retval = MPI_SUCCESS;
163   int src;
164   smpi_mpi_request_t *sendreq;
165   smpi_bench_end();
166   src = smpi_comm_rank(comm, SIMIX_host_self());
167   retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &sendreq);
168   if (NULL != sendreq) {
169     smpi_isend(sendreq);
170     smpi_wait(sendreq, MPI_STATUS_IGNORE);
171   }
172   smpi_bench_begin();
173   return retval;
174 }
175
176 int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm) {
177   int i;
178   int rank;
179   smpi_mpi_request_t *request;
180   smpi_bench_end();
181
182   rank = smpi_comm_rank(comm, SIMIX_host_self());
183
184   if (root == rank) {
185     smpi_create_request(buf, count, datatype, root, (rank + 1) % comm->size, 0, comm, &request);
186     if (comm->size > 2) {
187       request->fwdthrough = (rank - 1 + comm->size) % comm->size;
188     }
189     smpi_isend(request);
190     smpi_wait(request, MPI_STATUS_IGNORE);
191   } else {
192     smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank, 0, comm, &request);
193     if (NULL != request) {
194       smpi_irecv(request);
195       smpi_wait(request, MPI_STATUS_IGNORE);
196     }
197   }
198
199   smpi_bench_begin();
200   return MPI_SUCCESS;
201 }
202
203 int MPI_Barrier(MPI_Comm comm) {
204   smpi_bench_end();
205   smpi_barrier(comm);
206   smpi_bench_begin();
207   return MPI_SUCCESS;
208 }
209
210 // FIXME: instead of everyone sending in order, might be a good idea to send to next...
211 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
212   int i;
213   int rank;
214   smpi_mpi_request_t **sendreqs, **recvreqs;
215   smpi_bench_end();
216
217   rank = smpi_comm_rank(comm, SIMIX_host_self());
218
219   sendreqs = xbt_malloc(sizeof(smpi_mpi_request_t*) * comm->size);
220   recvreqs = xbt_malloc(sizeof(smpi_mpi_request_t*) * comm->size);
221
222   for (i = 0; i < comm->size; i++) {
223     if (rank == i) {
224       sendreqs[i] = NULL;
225       recvreqs[i] = NULL;
226       memcpy(recvbuf + recvtype->size * recvcount * i, sendbuf + sendtype->size * sendcount * i, recvtype->size * recvcount);
227     } else {
228       smpi_create_request(sendbuf + sendtype->size * sendcount * i, sendcount, sendtype, rank, i, 0, comm, sendreqs + i);
229       smpi_isend(sendreqs[i]);
230       smpi_create_request(recvbuf + recvtype->size * recvcount * i, recvcount, recvtype, i, rank, 0, comm, recvreqs + i);
231       smpi_irecv(recvreqs[i]);
232     }
233   } 
234
235   smpi_wait_all_nostatus(comm->size, sendreqs);
236   smpi_wait_all_nostatus(comm->size, recvreqs);
237
238   xbt_free(sendreqs);
239   xbt_free(recvreqs);
240
241   smpi_bench_begin();
242   return MPI_SUCCESS;
243 }
244
245 // FIXME: mpi routines shouldn't call mpi routines, complexity belongs at a lower level
246 // also, there's probably a really clever way to overlap everything for this one...
247 int MPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
248   MPI_Reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
249   MPI_Bcast(recvbuf, count, datatype, 0, comm);
250   return MPI_SUCCESS;
251 }
252
253 // FIXME: check if behavior defined when send/recv bufs are same...
254 int MPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) {
255   int i, j;
256   int rank;
257   smpi_mpi_request_t **requests;
258   void **scratchbuf;
259   smpi_bench_end();
260
261   rank = smpi_comm_rank(comm, SIMIX_host_self());
262
263   if (root == rank) {
264     requests = xbt_malloc(sizeof(smpi_mpi_request_t*) * comm->size);
265     scratchbuf = xbt_malloc(sizeof(void*) * comm->size);
266     memcpy(recvbuf, sendbuf, datatype->size * count);
267     for (i = 0; i < comm->size; i++) {
268       if (rank == i) {
269         requests[i] = NULL;
270         scratchbuf[i] = NULL;
271       } else {
272         scratchbuf[i] = xbt_malloc(datatype->size * count);
273         smpi_create_request(scratchbuf[i], count, datatype, MPI_ANY_SOURCE, rank, 0, comm, requests + i);
274         smpi_irecv(requests[i]);
275       }
276     }
277     smpi_wait_all_nostatus(comm->size, requests); // FIXME: use wait_any for slight performance gain
278     for (i = 0; i < comm->size; i++) {
279       if (rank != i) {
280         for (j = 0; j < count; j++) {
281           op->func(scratchbuf[i] + datatype->size * j, recvbuf + datatype->size * j, recvbuf + datatype->size * j);
282         }
283         xbt_free(requests[i]);
284         xbt_free(scratchbuf[i]);
285       }
286     }
287     xbt_free(requests);
288     xbt_free(scratchbuf);
289   } else {
290     requests = xbt_malloc(sizeof(smpi_mpi_request_t*));
291     smpi_create_request(sendbuf, count, datatype, rank, root, 0, comm, requests);
292     smpi_isend(*requests);
293     smpi_wait(*requests, MPI_STATUS_IGNORE);
294     xbt_free(*requests);
295     xbt_free(requests);
296   }
297
298   smpi_bench_begin();
299   return MPI_SUCCESS;
300 }