Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Implement MPI_Waitany and MPI_Waitall
[simgrid.git] / src / smpi / smpi_base.c
1 #include "private.h"
2 #include "xbt/time.h"
3
4 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_base, smpi,
5                                 "Logging specific to SMPI (base)");
6 XBT_LOG_EXTERNAL_CATEGORY(smpi_base);
7 XBT_LOG_EXTERNAL_CATEGORY(smpi_bench);
8 XBT_LOG_EXTERNAL_CATEGORY(smpi_kernel);
9 XBT_LOG_EXTERNAL_CATEGORY(smpi_mpi);
10 XBT_LOG_EXTERNAL_CATEGORY(smpi_receiver);
11 XBT_LOG_EXTERNAL_CATEGORY(smpi_sender);
12 XBT_LOG_EXTERNAL_CATEGORY(smpi_util);
13
14 smpi_mpi_global_t smpi_mpi_global = NULL;
15
16 /**
17  * Operations of MPI_OP : implemented=land,sum,min,max
18  **/
19 void smpi_mpi_land_func(void *a, void *b, int *length,
20                         MPI_Datatype * datatype);
21
22 void smpi_mpi_land_func(void *a, void *b, int *length,
23                         MPI_Datatype * datatype)
24 {
25   int i;
26   if (*datatype == smpi_mpi_global->mpi_int) {
27     int *x = a, *y = b;
28     for (i = 0; i < *length; i++) {
29       y[i] = x[i] && y[i];
30     }
31   }
32 }
33
34 /**
35  * sum two vectors element-wise
36  *
37  * @param a the first vectors
38  * @param b the second vectors
39  * @return the second vector is modified and contains the element-wise sums
40  **/
41 void smpi_mpi_sum_func(void *a, void *b, int *length,
42                        MPI_Datatype * datatype);
43
44 void smpi_mpi_sum_func(void *a, void *b, int *length, MPI_Datatype * datatype)
45 {
46           int i;
47           if (*datatype == smpi_mpi_global->mpi_byte) {
48                                 char *x = a, *y = b;
49                                 for (i = 0; i < *length; i++) {
50                                           y[i] = x[i] + y[i];
51                                 }
52           } else if (*datatype == smpi_mpi_global->mpi_int) {
53                                 int *x = a, *y = b;
54                                 for (i = 0; i < *length; i++) {
55                                           y[i] = x[i] + y[i];
56                                 }
57           } else if (*datatype == smpi_mpi_global->mpi_float) {
58                                 float *x = a, *y = b;
59                                 for (i = 0; i < *length; i++) {
60                                           y[i] = x[i] + y[i];
61                                 }
62           } else if (*datatype == smpi_mpi_global->mpi_double) {
63                                 double *x = a, *y = b;
64                                 for (i = 0; i < *length; i++) {
65                                           y[i] = x[i] + y[i];
66                                 }
67           }
68 }
69 /**
70  * compute the min of two vectors element-wise
71  **/
72 void smpi_mpi_min_func(void *a, void *b, int *length, MPI_Datatype * datatype);
73
74 void smpi_mpi_min_func(void *a, void *b, int *length, MPI_Datatype * datatype)
75 {
76           int i;
77           if (*datatype == smpi_mpi_global->mpi_byte) {
78                                 char *x = a, *y = b;
79                                 for (i = 0; i < *length; i++) {
80                                           y[i] = x[i] < y[i] ? x[i] : y[i];
81                                 }
82           } else {
83           if (*datatype == smpi_mpi_global->mpi_int) {
84                                 int *x = a, *y = b;
85                                 for (i = 0; i < *length; i++) {
86                                           y[i] = x[i] < y[i] ? x[i] : y[i];
87                                 }
88           } else {
89           if (*datatype == smpi_mpi_global->mpi_float) {
90                                 float *x = a, *y = b;
91                                 for (i = 0; i < *length; i++) {
92                                           y[i] = x[i] < y[i] ? x[i] : y[i];
93                                 }
94           } else {
95           if (*datatype == smpi_mpi_global->mpi_double) {
96                                 double *x = a, *y = b;
97                                 for (i = 0; i < *length; i++) {
98                                           y[i] = x[i] < y[i] ? x[i] : y[i];
99                                 }
100
101           }}}}
102 }
103 /**
104  * compute the max of two vectors element-wise
105  **/
106 void smpi_mpi_max_func(void *a, void *b, int *length, MPI_Datatype * datatype);
107
108 void smpi_mpi_max_func(void *a, void *b, int *length, MPI_Datatype * datatype)
109 {
110           int i;
111           if (*datatype == smpi_mpi_global->mpi_byte) {
112                                 char *x = a, *y = b;
113                                 for (i = 0; i < *length; i++) {
114                                           y[i] = x[i] > y[i] ? x[i] : y[i];
115                                 }
116           } else {
117           if (*datatype == smpi_mpi_global->mpi_int) {
118                                 int *x = a, *y = b;
119                                 for (i = 0; i > *length; i++) {
120                                           y[i] = x[i] < y[i] ? x[i] : y[i];
121                                 }
122           } else {
123           if (*datatype == smpi_mpi_global->mpi_float) {
124                                 float *x = a, *y = b;
125                                 for (i = 0; i > *length; i++) {
126                                           y[i] = x[i] < y[i] ? x[i] : y[i];
127                                 }
128           } else {
129           if (*datatype == smpi_mpi_global->mpi_double) {
130                                 double *x = a, *y = b;
131                                 for (i = 0; i > *length; i++) {
132                                           y[i] = x[i] < y[i] ? x[i] : y[i];
133                                 }
134
135           }}}}
136 }
137
138
139
140
141 /**
142  * tell the MPI rank of the calling process (from its SIMIX process id)
143  **/
144 int smpi_mpi_comm_rank(smpi_mpi_communicator_t comm)
145 {
146   return comm->index_to_rank_map[smpi_process_index()];
147 }
148
149 void smpi_process_init(int *argc, char***argv)
150 {
151   smpi_process_data_t pdata;
152
153   // initialize some local variables
154
155   pdata = xbt_new(s_smpi_process_data_t, 1);
156   SIMIX_process_set_data(SIMIX_process_self(),pdata);
157
158   /* get rank from command line, and remove it from argv */
159   pdata->index = atoi( (*argv)[1] );
160   DEBUG1("I'm rank %d",pdata->index);
161   if (*argc>2) {
162           memmove((*argv)[1],(*argv)[2], sizeof(char*)* (*argc-2));
163           (*argv)[ (*argc)-1] = NULL;
164   }
165   (*argc)--;
166
167   pdata->mutex = SIMIX_mutex_init();
168   pdata->cond = SIMIX_cond_init();
169   pdata->finalize = 0;
170
171   pdata->pending_recv_request_queue = xbt_fifo_new();
172   pdata->pending_send_request_queue = xbt_fifo_new();
173   pdata->received_message_queue = xbt_fifo_new();
174
175   pdata->main = SIMIX_process_self();
176   pdata->sender = SIMIX_process_create("smpi_sender",
177           smpi_sender, pdata,
178           SIMIX_host_get_name(SIMIX_host_self()), 0, NULL,
179           /*props */ NULL);
180   pdata->receiver = SIMIX_process_create("smpi_receiver",
181           smpi_receiver, pdata,
182           SIMIX_host_get_name(SIMIX_host_self()), 0, NULL,
183           /*props */ NULL);
184
185   smpi_global->main_processes[pdata->index] = SIMIX_process_self();
186   return;
187 }
188
189 void smpi_process_finalize()
190 {
191   smpi_process_data_t pdata =  SIMIX_process_get_data(SIMIX_process_self());
192
193   pdata->finalize = 2; /* Tell sender and receiver to quit */
194   SIMIX_process_resume(pdata->sender);
195   SIMIX_process_resume(pdata->receiver);
196   while (pdata->finalize>0) { /* wait until it's done */
197           SIMIX_cond_wait(pdata->cond,pdata->mutex);
198   }
199
200   SIMIX_mutex_destroy(pdata->mutex);
201   SIMIX_cond_destroy(pdata->cond);
202   xbt_fifo_free(pdata->pending_recv_request_queue);
203   xbt_fifo_free(pdata->pending_send_request_queue);
204   xbt_fifo_free(pdata->received_message_queue);
205 }
206
207 int smpi_mpi_barrier(smpi_mpi_communicator_t comm)
208 {
209
210   SIMIX_mutex_lock(comm->barrier_mutex);
211   ++comm->barrier_count;
212   if (comm->barrier_count > comm->size) {       // only happens on second barrier...
213     comm->barrier_count = 0;
214   } else if (comm->barrier_count == comm->size) {
215     SIMIX_cond_broadcast(comm->barrier_cond);
216   }
217   while (comm->barrier_count < comm->size) {
218     SIMIX_cond_wait(comm->barrier_cond, comm->barrier_mutex);
219   }
220   SIMIX_mutex_unlock(comm->barrier_mutex);
221
222   return MPI_SUCCESS;
223 }
224
225 int smpi_mpi_isend(smpi_mpi_request_t request)
226 {
227         smpi_process_data_t pdata =  SIMIX_process_get_data(SIMIX_process_self());
228   int retval = MPI_SUCCESS;
229
230   if (NULL == request) {
231     retval = MPI_ERR_INTERN;
232   } else {
233     xbt_fifo_push(pdata->pending_send_request_queue, request);
234     SIMIX_process_resume(pdata->sender);
235   }
236
237   return retval;
238 }
239
240 int smpi_mpi_irecv(smpi_mpi_request_t request)
241 {
242   int retval = MPI_SUCCESS;
243   smpi_process_data_t pdata =  SIMIX_process_get_data(SIMIX_process_self());
244
245   if (NULL == request) {
246     retval = MPI_ERR_INTERN;
247   } else {
248     xbt_fifo_push(pdata->pending_recv_request_queue, request);
249
250     if (SIMIX_process_is_suspended(pdata->receiver)) {
251       SIMIX_process_resume(pdata->receiver);
252     }
253   }
254
255   return retval;
256 }
257
258 int smpi_mpi_wait(smpi_mpi_request_t request, smpi_mpi_status_t * status)
259 {
260   int retval = MPI_SUCCESS;
261
262   if (NULL == request) {
263     retval = MPI_ERR_INTERN;
264   } else {
265     SIMIX_mutex_lock(request->mutex);
266     while (!request->completed) {
267       SIMIX_cond_wait(request->cond, request->mutex);
268     }
269     if (NULL != status) {
270       status->MPI_SOURCE = request->src;
271       status->MPI_TAG = request->tag;
272       status->MPI_ERROR = MPI_SUCCESS;
273     }
274     SIMIX_mutex_unlock(request->mutex);
275   }
276
277   return retval;
278 }
279
280 int smpi_mpi_wait_all(int count, smpi_mpi_request_t *requests, smpi_mpi_status_t **status) {
281         int cpt;
282         int index;
283         int retval;
284         smpi_mpi_status_t stat;
285
286         for (cpt=0; cpt<count;cpt++) {
287                 retval = smpi_mpi_wait_any(count,requests, &index,&stat);
288                 if (retval != MPI_SUCCESS)
289                         return retval;
290                 memcpy(status[index],&stat,sizeof(stat));
291         }
292         return MPI_SUCCESS;
293 }
294
295 int smpi_mpi_wait_any(int count, smpi_mpi_request_t *requests, int *index, smpi_mpi_status_t *status) {
296           int cpt;
297
298           *index = MPI_UNDEFINED;
299           if (NULL == requests) {
300             return MPI_ERR_INTERN;
301           }
302           /* First check if one of them is already done */
303           for (cpt=0;cpt<count;cpt++) {
304                   if (requests[cpt]->completed && !requests[cpt]->consumed) { /* got ya */
305                           *index=cpt;
306                           goto found_request;
307                   }
308           }
309           /* If none found, block */
310           /* FIXME: should use a SIMIX_cond_waitany, when implemented. For now, block on the first one */
311           while (1) {
312                   for (cpt=0;cpt<count;cpt++) {
313                           if (!requests[cpt]->completed) { /* this one is not done, wait on it */
314                                   while (!requests[cpt]->completed)
315                                       SIMIX_cond_wait(requests[cpt]->cond, requests[cpt]->mutex);
316
317                                   *index=cpt;
318                                   goto found_request;
319                           }
320                   }
321                   if (cpt == count) /* they are all done. Damn user */
322                           return MPI_ERR_REQUEST;
323           }
324
325           found_request:
326           requests[*index]->consumed = 1;
327
328            if (NULL != status) {
329               status->MPI_SOURCE = requests[*index]->src;
330               status->MPI_TAG = requests[*index]->tag;
331               status->MPI_ERROR = MPI_SUCCESS;
332             }
333           return MPI_SUCCESS;
334
335 }