Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
4ebd97af96d61bb061d945edcf7e13deaffcfd01
[simgrid.git] / src / smpi / smpi_global.c
1 /* Copyright (c) 2007-2013. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5   * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include <stdint.h>
8 #include <stdio.h>
9 #include <stdlib.h>
10
11 #include "private.h"
12 #include "smpi_mpi_dt_private.h"
13 #include "mc/mc.h"
14 #include "surf/surf.h"
15 #include "simix/smx_private.h"
16 #include "simgrid/sg_config.h"
17
18
19 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_kernel, smpi,
20                                 "Logging specific to SMPI (kernel)");
21
22 typedef struct s_smpi_process_data {
23   int index;
24   int* argc;
25   char*** argv;
26   smx_rdv_t mailbox;
27   smx_rdv_t mailbox_small;
28   xbt_os_timer_t timer;
29   double simulated;
30   MPI_Comm comm_self;
31   void *data; /* user data */
32   int initialized;
33 } s_smpi_process_data_t;
34
35 static smpi_process_data_t *process_data = NULL;
36 static int process_count = 0;
37
38 MPI_Comm MPI_COMM_WORLD = MPI_COMM_NULL;
39 int MPI_UNIVERSE_SIZE;
40
41 MPI_Errhandler* MPI_ERRORS_RETURN = NULL;
42 MPI_Errhandler* MPI_ERRORS_ARE_FATAL = NULL;
43 MPI_Errhandler* MPI_ERRHANDLER_NULL = NULL;
44
45 #define MAILBOX_NAME_MAXLEN (5 + sizeof(int) * 2 + 1)
46
47 static char* get_mailbox_name(char* str, int index) {
48   snprintf(str, MAILBOX_NAME_MAXLEN, "SMPI-%0*x", (int)(sizeof(int) * 2), index);
49   return str;
50 }
51
52 static char* get_mailbox_name_small(char* str, int index) {
53   snprintf(str, MAILBOX_NAME_MAXLEN, "small%0*x", (int)(sizeof(int) * 2), index);
54   return str;
55 }
56
57 void smpi_process_init(int *argc, char ***argv)
58 {
59   int index;
60   smpi_process_data_t data;
61   smx_process_t proc;
62
63   if(argc && argv) {
64     proc = SIMIX_process_self();
65     index = atoi((*argv)[1]);
66     data = smpi_process_remote_data(index);
67     simcall_process_set_data(proc, data);
68     if (*argc > 2) {
69       free((*argv)[1]);
70       memmove(&(*argv)[1], &(*argv)[2], sizeof(char *) * (*argc - 2));
71       (*argv)[(*argc) - 1] = NULL;
72     }
73     (*argc)--;
74     data->argc = argc;
75     data->argv = argv;
76     simcall_rdv_set_receiver(data->mailbox_small, proc);// set the process attached to the mailbox
77     XBT_DEBUG("<%d> New process in the game: %p", index, proc);
78   }
79 }
80
81 void smpi_process_destroy(void)
82 {
83   int index = smpi_process_index();
84   process_data[index]->index=-100;
85   XBT_DEBUG("<%d> Process left the game", index);
86 }
87
88 /**
89  * @brief Prepares the current process for termination.
90  */
91 void smpi_process_finalize(void)
92 {
93   // wait for all pending asynchronous comms to finish
94   while (SIMIX_process_has_pending_comms(SIMIX_process_self())) {
95     simcall_process_sleep(0.01);
96   }
97 }
98
99 /**
100  * @brief Check if a process is finalized
101  */
102 int smpi_process_finalized()
103 {
104    return (smpi_process_index()==-100);
105   // If finalized, this value has been set to -100;
106 }
107
108 /**
109  * @brief Check if a process is initialized
110  */
111 int smpi_process_initialized(void)
112 {
113   int index = smpi_process_index();
114   return((index != -100) && (index!=MPI_UNDEFINED) && (process_data[index]->initialized));
115 }
116
117 /**
118  * @brief Mark a process as initialized (=MPI_Init called)
119  */
120 void smpi_process_mark_as_initialized(void)
121 {
122   int index = smpi_process_index();
123   if((index != -100)&& (index!=MPI_UNDEFINED))process_data[index]->initialized=1;
124 }
125
126
127 #ifdef SMPI_F2C
128 int smpi_process_argc(void) {
129   smpi_process_data_t data = smpi_process_data();
130
131   return data->argc ? *(data->argc) - 1 : 0;
132 }
133
134 int smpi_process_getarg(integer* index, char* dst, ftnlen len) {
135   smpi_process_data_t data = smpi_process_data();
136   char* arg;
137   ftnlen i;
138
139   if(!data->argc || !data->argv
140      || *index < 1 || *index >= *(data->argc)) {
141     return -1;
142   }
143   arg = (*data->argv)[*index];
144   for(i = 0; i < len && arg[i] != '\0'; i++) {
145     dst[i] = arg[i];
146   }
147   for(; i < len; i++) {
148     dst[i] = ' ';
149   }
150   return 0;
151 }
152
153 int smpi_global_size(void) {
154    char* value = getenv("SMPI_GLOBAL_SIZE");
155
156    if(!value) {
157      fprintf(stderr, "Please set env var SMPI_GLOBAL_SIZE to expected number of processes.\n");
158      xbt_abort();
159    }
160    return atoi(value);
161 }
162 #endif
163
164 smpi_process_data_t smpi_process_data(void)
165 {
166   return SIMIX_process_self_get_data(SIMIX_process_self());
167 }
168
169 smpi_process_data_t smpi_process_remote_data(int index)
170 {
171   return process_data[index];
172 }
173
174 void smpi_process_set_user_data(void *data)
175 {
176   smpi_process_data_t process_data = smpi_process_data();
177   process_data->data = data;
178 }
179
180 void* smpi_process_get_user_data(){
181   smpi_process_data_t process_data = smpi_process_data();
182   return process_data->data;
183 }
184
185 int smpi_process_count(void)
186 {
187   return process_count;
188 }
189
190 int smpi_process_index(void)
191 {
192   smpi_process_data_t data = smpi_process_data();
193   //return -1 if not initialized
194   return data? data->index : MPI_UNDEFINED;
195 }
196
197 smx_rdv_t smpi_process_mailbox(void) {
198   smpi_process_data_t data = smpi_process_data();
199
200   return data->mailbox;
201 }
202
203 smx_rdv_t smpi_process_mailbox_small(void) {
204   smpi_process_data_t data = smpi_process_data();
205
206   return data->mailbox_small;
207 }
208
209 smx_rdv_t smpi_process_remote_mailbox(int index) {
210   smpi_process_data_t data = smpi_process_remote_data(index);
211
212   return data->mailbox;
213 }
214
215
216 smx_rdv_t smpi_process_remote_mailbox_small(int index) {
217   smpi_process_data_t data = smpi_process_remote_data(index);
218
219   return data->mailbox_small;
220 }
221
222 xbt_os_timer_t smpi_process_timer(void)
223 {
224   smpi_process_data_t data = smpi_process_data();
225
226   return data->timer;
227 }
228
229 void smpi_process_simulated_start(void) {
230   smpi_process_data_t data = smpi_process_data();
231
232   data->simulated = SIMIX_get_clock();
233 }
234
235 double smpi_process_simulated_elapsed(void) {
236   smpi_process_data_t data = smpi_process_data();
237
238   return SIMIX_get_clock() - data->simulated;
239 }
240
241 MPI_Comm smpi_process_comm_self(void) {
242   smpi_process_data_t data = smpi_process_data();
243
244   return data->comm_self;
245 }
246
247 void print_request(const char *message, MPI_Request request) {
248   XBT_DEBUG("%s  request %p  [buf = %p, size = %zu, src = %d, dst = %d, tag = %d, flags = %x]",
249          message, request, request->buf, request->size,
250          request->src, request->dst, request->tag, request->flags);
251 }
252
253 static void SMPI_comm_copy_buffer_callback(smx_action_t comm, void* buff, size_t buff_size)
254 {
255   XBT_DEBUG("Copy the data over");
256   memcpy(comm->comm.dst_buff, buff, buff_size);
257   if (comm->comm.detached) { // if this is a detached send, the source buffer was duplicated by SMPI sender to make the original buffer available to the application ASAP
258     xbt_free(buff);
259     //It seems that the request is used after the call there this should
260     //be free somewhereelse  but where???
261     //xbt_free(comm->comm.src_data);// inside SMPI the request is keep
262     //inside the user data and should be free 
263     comm->comm.src_buff = NULL;
264   }
265 }
266
267 void smpi_global_init(void)
268 {
269   int i;
270   MPI_Group group;
271   char name[MAILBOX_NAME_MAXLEN];
272
273   SIMIX_comm_set_copy_data_callback(&SMPI_comm_copy_buffer_callback);
274   process_count = SIMIX_process_count();
275   process_data = xbt_new(smpi_process_data_t, process_count);
276   for (i = 0; i < process_count; i++) {
277     process_data[i] = xbt_new(s_smpi_process_data_t, 1);
278     process_data[i]->index = i;
279     process_data[i]->argc = NULL;
280     process_data[i]->argv = NULL;
281     process_data[i]->mailbox = simcall_rdv_create(get_mailbox_name(name, i));
282     process_data[i]->mailbox_small = simcall_rdv_create(get_mailbox_name_small(name, i));
283     process_data[i]->timer = xbt_os_timer_new();
284     if(MC_is_active())
285       MC_ignore_heap(process_data[i]->timer, xbt_os_timer_size());
286     group = smpi_group_new(1);
287     process_data[i]->comm_self = smpi_comm_new(group);
288     process_data[i]->initialized =0;
289
290     smpi_group_set_mapping(group, i, 0);
291   }
292   group = smpi_group_new(process_count);
293   MPI_COMM_WORLD = smpi_comm_new(group);
294   MPI_UNIVERSE_SIZE = smpi_comm_size(MPI_COMM_WORLD);
295   for (i = 0; i < process_count; i++) {
296     smpi_group_set_mapping(group, i, i);
297   }
298
299   //check correctness of MPI parameters
300
301   xbt_assert(sg_cfg_get_int("smpi/async_small_thres")<=sg_cfg_get_int("smpi/send_is_detached_thres"));
302 }
303
304 void smpi_global_destroy(void)
305 {
306   int count = smpi_process_count();
307   int i;
308
309   smpi_bench_destroy();
310   while(smpi_group_unuse(smpi_comm_group(MPI_COMM_WORLD))>0);
311   xbt_free(MPI_COMM_WORLD);
312   MPI_COMM_WORLD = MPI_COMM_NULL;
313   for (i = 0; i < count; i++) {
314     smpi_group_unuse(smpi_comm_group(process_data[i]->comm_self));
315     smpi_comm_destroy(process_data[i]->comm_self);
316     xbt_os_timer_free(process_data[i]->timer);
317     simcall_rdv_destroy(process_data[i]->mailbox);
318     simcall_rdv_destroy(process_data[i]->mailbox_small);
319     xbt_free(process_data[i]);
320   }
321   xbt_free(process_data);
322   process_data = NULL;
323
324   smpi_free_static();
325 }
326
327 /* Fortran specific stuff */
328 /* With smpicc, the following weak symbols are used */
329 /* With smpiff, the following weak symbols are replaced by those in libf2c */
330 int __attribute__((weak)) xargc;
331 char** __attribute__((weak)) xargv;
332
333 #ifndef WIN32
334 void __attribute__((weak)) user_main_(){
335   xbt_die("Should not be in this smpi_simulated_main");
336   return;
337 }
338 int __attribute__((weak)) smpi_simulated_main_(int argc, char** argv) {
339   smpi_process_init(&argc, &argv);
340   user_main_();
341   //xbt_die("Should not be in this smpi_simulated_main");
342   return 0;
343 }
344
345 int __attribute__((weak)) main(int argc, char** argv) {
346    return smpi_main(smpi_simulated_main_,argc,argv);
347 }
348
349 int __attribute__((weak)) MAIN__(){
350   return smpi_main(smpi_simulated_main_,xargc, xargv);
351 };
352 #endif
353
354 int smpi_main(int (*realmain) (int argc, char *argv[]),int argc, char *argv[])
355 {
356   srand(SMPI_RAND_SEED);
357   
358   if(getenv("SMPI_PRETEND_CC") != NULL) {
359   /* Hack to ensure that smpicc can pretend to be a simple compiler. Particularly handy to pass it to the configuration tools */
360     return 0;
361   }
362
363   /* Connect log categories.  See xbt/log.c */
364   XBT_LOG_CONNECT(smpi);  /* Keep this line as soon as possible in this
365                              function: xbt_log_appender_file.c depends on it
366                              DO NOT connect this in XBT or so, or it will be
367                              useless to xbt_log_appender_file.c */
368 #ifdef HAVE_TRACING
369   XBT_LOG_CONNECT(instr_smpi);
370 #endif
371   XBT_LOG_CONNECT(smpi_base);
372   XBT_LOG_CONNECT(smpi_bench);
373   XBT_LOG_CONNECT(smpi_coll);
374   XBT_LOG_CONNECT(smpi_comm);
375   XBT_LOG_CONNECT(smpi_group);
376   XBT_LOG_CONNECT(smpi_kernel);
377   XBT_LOG_CONNECT(smpi_mpi);
378   XBT_LOG_CONNECT(smpi_mpi_dt);
379   XBT_LOG_CONNECT(smpi_pmpi);
380   XBT_LOG_CONNECT(smpi_replay);
381   XBT_LOG_CONNECT(smpi_colls);
382
383 #ifdef HAVE_TRACING
384   TRACE_global_init(&argc, argv);
385
386   TRACE_add_start_function(TRACE_smpi_alloc);
387   TRACE_add_end_function(TRACE_smpi_release);
388 #endif
389
390   SIMIX_global_init(&argc, argv);
391
392 #ifdef HAVE_TRACING
393   TRACE_start();
394 #endif
395
396   // parse the platform file: get the host list
397   SIMIX_create_environment(argv[1]);
398
399   SIMIX_function_register_default(realmain);
400   SIMIX_launch_application(argv[2]);
401
402   int gather_id = find_coll_description(mpi_coll_gather_description,
403                                            sg_cfg_get_string("smpi/gather"));
404   mpi_coll_gather_fun = (int (*)(void *, int, MPI_Datatype,
405                                     void*, int, MPI_Datatype, int, MPI_Comm))
406                            mpi_coll_gather_description[gather_id].coll;
407
408   int allgather_id = find_coll_description(mpi_coll_allgather_description,
409                                            sg_cfg_get_string("smpi/allgather"));
410   mpi_coll_allgather_fun = (int (*)(void *, int, MPI_Datatype,
411                                     void*, int, MPI_Datatype, MPI_Comm))
412                            mpi_coll_allgather_description[allgather_id].coll;
413
414   int allgatherv_id = find_coll_description(mpi_coll_allgatherv_description,
415                                            sg_cfg_get_string("smpi/allgatherv"));
416   mpi_coll_allgatherv_fun = (int (*)(void *, int, MPI_Datatype,
417                                     void*, int*, int*, MPI_Datatype, MPI_Comm))
418                            mpi_coll_allgatherv_description[allgatherv_id].coll;
419
420   int allreduce_id = find_coll_description(mpi_coll_allreduce_description,
421                                            sg_cfg_get_string("smpi/allreduce"));
422   mpi_coll_allreduce_fun = (int (*)(void *sbuf, void *rbuf, int rcount, \
423                                     MPI_Datatype dtype, MPI_Op op, MPI_Comm comm))
424                            mpi_coll_allreduce_description[allreduce_id].coll;
425
426   int alltoall_id = find_coll_description(mpi_coll_alltoall_description,
427                                           sg_cfg_get_string("smpi/alltoall"));
428   mpi_coll_alltoall_fun = (int (*)(void *, int, MPI_Datatype,
429                                    void*, int, MPI_Datatype, MPI_Comm))
430                           mpi_coll_alltoall_description[alltoall_id].coll;
431
432   int alltoallv_id = find_coll_description(mpi_coll_alltoallv_description,
433                                           sg_cfg_get_string("smpi/alltoallv"));
434   mpi_coll_alltoallv_fun = (int (*)(void *, int*, int*, MPI_Datatype,
435                                     void*, int*, int*, MPI_Datatype, MPI_Comm))
436                           mpi_coll_alltoallv_description[alltoallv_id].coll;
437
438   int bcast_id = find_coll_description(mpi_coll_bcast_description,
439                                           sg_cfg_get_string("smpi/bcast"));
440   mpi_coll_bcast_fun = (int (*)(void *buf, int count, MPI_Datatype datatype, \
441                                 int root, MPI_Comm com))
442                        mpi_coll_bcast_description[bcast_id].coll;
443
444   int reduce_id = find_coll_description(mpi_coll_reduce_description,
445                                           sg_cfg_get_string("smpi/reduce"));
446   mpi_coll_reduce_fun = (int (*)(void *buf, void *rbuf, int count, MPI_Datatype datatype, \
447                                  MPI_Op op, int root, MPI_Comm comm))
448                         mpi_coll_reduce_description[reduce_id].coll;
449
450   int reduce_scatter_id = find_coll_description(mpi_coll_reduce_scatter_description,
451                                            sg_cfg_get_string("smpi/reduce_scatter"));
452   mpi_coll_reduce_scatter_fun = (int (*)(void *sbuf, void *rbuf, int *rcounts,\
453                     MPI_Datatype dtype,MPI_Op  op,MPI_Comm  comm))
454                            mpi_coll_reduce_scatter_description[reduce_scatter_id].coll;
455
456   int scatter_id = find_coll_description(mpi_coll_scatter_description,
457                                            sg_cfg_get_string("smpi/scatter"));
458   mpi_coll_scatter_fun = (int (*)(void *sendbuf, int sendcount, MPI_Datatype sendtype,\
459                 void *recvbuf, int recvcount, MPI_Datatype recvtype,\
460                 int root, MPI_Comm comm))
461                            mpi_coll_scatter_description[scatter_id].coll;
462
463   int barrier_id = find_coll_description(mpi_coll_barrier_description,
464                                            sg_cfg_get_string("smpi/barrier"));
465   mpi_coll_barrier_fun = (int (*)(MPI_Comm comm))
466                            mpi_coll_barrier_description[barrier_id].coll;
467
468   smpi_global_init();
469
470   /* Clean IO before the run */
471   fflush(stdout);
472   fflush(stderr);
473
474   if (MC_is_active())
475     MC_do_the_modelcheck_for_real();
476   else
477     SIMIX_run();
478
479   if (sg_cfg_get_boolean("smpi/display_timing"))
480     XBT_INFO("Simulation time: %g seconds.", SIMIX_get_clock());
481
482   smpi_global_destroy();
483
484 #ifdef HAVE_TRACING
485   TRACE_end();
486 #endif
487
488   return 0;
489 }